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 #include <petsc/private/vecimpl.h> 13 #include <petsc/private/hashseti.h> 14 15 #define MATIS_MAX_ENTRIES_INSERTION 2048 16 static PetscErrorCode MatSetValuesLocal_IS(Mat,PetscInt,const PetscInt*,PetscInt,const PetscInt*,const PetscScalar*,InsertMode); 17 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat,PetscInt,const PetscInt*,PetscInt,const PetscInt*,const PetscScalar*,InsertMode); 18 static PetscErrorCode MatISSetUpScatters_Private(Mat); 19 20 static PetscErrorCode MatISContainerDestroyPtAP_Private(void *ptr) 21 { 22 MatISPtAP ptap = (MatISPtAP)ptr; 23 24 PetscFunctionBegin; 25 PetscCall(MatDestroySubMatrices(ptap->ris1 ? 2 : 1,&ptap->lP)); 26 PetscCall(ISDestroy(&ptap->cis0)); 27 PetscCall(ISDestroy(&ptap->cis1)); 28 PetscCall(ISDestroy(&ptap->ris0)); 29 PetscCall(ISDestroy(&ptap->ris1)); 30 PetscCall(PetscFree(ptap)); 31 PetscFunctionReturn(0); 32 } 33 34 static PetscErrorCode MatPtAPNumeric_IS_XAIJ(Mat A, Mat P, Mat C) 35 { 36 MatISPtAP ptap; 37 Mat_IS *matis = (Mat_IS*)A->data; 38 Mat lA,lC; 39 MatReuse reuse; 40 IS ris[2],cis[2]; 41 PetscContainer c; 42 PetscInt n; 43 44 PetscFunctionBegin; 45 PetscCall(PetscObjectQuery((PetscObject)C,"_MatIS_PtAP",(PetscObject*)&c)); 46 PetscCheck(c,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing PtAP information"); 47 PetscCall(PetscContainerGetPointer(c,(void**)&ptap)); 48 ris[0] = ptap->ris0; 49 ris[1] = ptap->ris1; 50 cis[0] = ptap->cis0; 51 cis[1] = ptap->cis1; 52 n = ptap->ris1 ? 2 : 1; 53 reuse = ptap->lP ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 54 PetscCall(MatCreateSubMatrices(P,n,ris,cis,reuse,&ptap->lP)); 55 56 PetscCall(MatISGetLocalMat(A,&lA)); 57 PetscCall(MatISGetLocalMat(C,&lC)); 58 if (ptap->ris1) { /* unsymmetric A mapping */ 59 Mat lPt; 60 61 PetscCall(MatTranspose(ptap->lP[1],MAT_INITIAL_MATRIX,&lPt)); 62 PetscCall(MatMatMatMult(lPt,lA,ptap->lP[0],reuse,ptap->fill,&lC)); 63 if (matis->storel2l) { 64 PetscCall(PetscObjectCompose((PetscObject)(A),"_MatIS_PtAP_l2l",(PetscObject)lPt)); 65 } 66 PetscCall(MatDestroy(&lPt)); 67 } else { 68 PetscCall(MatPtAP(lA,ptap->lP[0],reuse,ptap->fill,&lC)); 69 if (matis->storel2l) { 70 PetscCall(PetscObjectCompose((PetscObject)C,"_MatIS_PtAP_l2l",(PetscObject)ptap->lP[0])); 71 } 72 } 73 if (reuse == MAT_INITIAL_MATRIX) { 74 PetscCall(MatISSetLocalMat(C,lC)); 75 PetscCall(MatDestroy(&lC)); 76 } 77 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 78 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 79 PetscFunctionReturn(0); 80 } 81 82 static PetscErrorCode MatGetNonzeroColumnsLocal_Private(Mat PT,IS *cis) 83 { 84 Mat Po,Pd; 85 IS zd,zo; 86 const PetscInt *garray; 87 PetscInt *aux,i,bs; 88 PetscInt dc,stc,oc,ctd,cto; 89 PetscBool ismpiaij,ismpibaij,isseqaij,isseqbaij; 90 MPI_Comm comm; 91 92 PetscFunctionBegin; 93 PetscValidHeaderSpecific(PT,MAT_CLASSID,1); 94 PetscValidPointer(cis,2); 95 PetscCall(PetscObjectGetComm((PetscObject)PT,&comm)); 96 bs = 1; 97 PetscCall(PetscObjectBaseTypeCompare((PetscObject)PT,MATMPIAIJ,&ismpiaij)); 98 PetscCall(PetscObjectBaseTypeCompare((PetscObject)PT,MATMPIBAIJ,&ismpibaij)); 99 PetscCall(PetscObjectBaseTypeCompare((PetscObject)PT,MATSEQAIJ,&isseqaij)); 100 PetscCall(PetscObjectTypeCompare((PetscObject)PT,MATSEQBAIJ,&isseqbaij)); 101 if (isseqaij || isseqbaij) { 102 Pd = PT; 103 Po = NULL; 104 garray = NULL; 105 } else if (ismpiaij) { 106 PetscCall(MatMPIAIJGetSeqAIJ(PT,&Pd,&Po,&garray)); 107 } else if (ismpibaij) { 108 PetscCall(MatMPIBAIJGetSeqBAIJ(PT,&Pd,&Po,&garray)); 109 PetscCall(MatGetBlockSize(PT,&bs)); 110 } else SETERRQ(comm,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(PT))->type_name); 111 112 /* identify any null columns in Pd or Po */ 113 /* We use a tolerance comparison since it may happen that, with geometric multigrid, 114 some of the columns are not really zero, but very close to */ 115 zo = zd = NULL; 116 if (Po) { 117 PetscCall(MatFindNonzeroRowsOrCols_Basic(Po,PETSC_TRUE,PETSC_SMALL,&zo)); 118 } 119 PetscCall(MatFindNonzeroRowsOrCols_Basic(Pd,PETSC_TRUE,PETSC_SMALL,&zd)); 120 121 PetscCall(MatGetLocalSize(PT,NULL,&dc)); 122 PetscCall(MatGetOwnershipRangeColumn(PT,&stc,NULL)); 123 if (Po) PetscCall(MatGetLocalSize(Po,NULL,&oc)); 124 else oc = 0; 125 PetscCall(PetscMalloc1((dc+oc)/bs,&aux)); 126 if (zd) { 127 const PetscInt *idxs; 128 PetscInt nz; 129 130 /* this will throw an error if bs is not valid */ 131 PetscCall(ISSetBlockSize(zd,bs)); 132 PetscCall(ISGetLocalSize(zd,&nz)); 133 PetscCall(ISGetIndices(zd,&idxs)); 134 ctd = nz/bs; 135 for (i=0; i<ctd; i++) aux[i] = (idxs[bs*i]+stc)/bs; 136 PetscCall(ISRestoreIndices(zd,&idxs)); 137 } else { 138 ctd = dc/bs; 139 for (i=0; i<ctd; i++) aux[i] = i+stc/bs; 140 } 141 if (zo) { 142 const PetscInt *idxs; 143 PetscInt nz; 144 145 /* this will throw an error if bs is not valid */ 146 PetscCall(ISSetBlockSize(zo,bs)); 147 PetscCall(ISGetLocalSize(zo,&nz)); 148 PetscCall(ISGetIndices(zo,&idxs)); 149 cto = nz/bs; 150 for (i=0; i<cto; i++) aux[i+ctd] = garray[idxs[bs*i]/bs]; 151 PetscCall(ISRestoreIndices(zo,&idxs)); 152 } else { 153 cto = oc/bs; 154 for (i=0; i<cto; i++) aux[i+ctd] = garray[i]; 155 } 156 PetscCall(ISCreateBlock(comm,bs,ctd+cto,aux,PETSC_OWN_POINTER,cis)); 157 PetscCall(ISDestroy(&zd)); 158 PetscCall(ISDestroy(&zo)); 159 PetscFunctionReturn(0); 160 } 161 162 static PetscErrorCode MatPtAPSymbolic_IS_XAIJ(Mat A,Mat P,PetscReal fill,Mat C) 163 { 164 Mat PT,lA; 165 MatISPtAP ptap; 166 ISLocalToGlobalMapping Crl2g,Ccl2g,rl2g,cl2g; 167 PetscContainer c; 168 MatType lmtype; 169 const PetscInt *garray; 170 PetscInt ibs,N,dc; 171 MPI_Comm comm; 172 173 PetscFunctionBegin; 174 PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 175 PetscCall(MatSetType(C,MATIS)); 176 PetscCall(MatISGetLocalMat(A,&lA)); 177 PetscCall(MatGetType(lA,&lmtype)); 178 PetscCall(MatISSetLocalMatType(C,lmtype)); 179 PetscCall(MatGetSize(P,NULL,&N)); 180 PetscCall(MatGetLocalSize(P,NULL,&dc)); 181 PetscCall(MatSetSizes(C,dc,dc,N,N)); 182 /* Not sure about this 183 PetscCall(MatGetBlockSizes(P,NULL,&ibs)); 184 PetscCall(MatSetBlockSize(*C,ibs)); 185 */ 186 187 PetscCall(PetscNew(&ptap)); 188 PetscCall(PetscContainerCreate(PETSC_COMM_SELF,&c)); 189 PetscCall(PetscContainerSetPointer(c,ptap)); 190 PetscCall(PetscContainerSetUserDestroy(c,MatISContainerDestroyPtAP_Private)); 191 PetscCall(PetscObjectCompose((PetscObject)C,"_MatIS_PtAP",(PetscObject)c)); 192 PetscCall(PetscContainerDestroy(&c)); 193 ptap->fill = fill; 194 195 PetscCall(MatISGetLocalToGlobalMapping(A,&rl2g,&cl2g)); 196 197 PetscCall(ISLocalToGlobalMappingGetBlockSize(cl2g,&ibs)); 198 PetscCall(ISLocalToGlobalMappingGetSize(cl2g,&N)); 199 PetscCall(ISLocalToGlobalMappingGetBlockIndices(cl2g,&garray)); 200 PetscCall(ISCreateBlock(comm,ibs,N/ibs,garray,PETSC_COPY_VALUES,&ptap->ris0)); 201 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(cl2g,&garray)); 202 203 PetscCall(MatCreateSubMatrix(P,ptap->ris0,NULL,MAT_INITIAL_MATRIX,&PT)); 204 PetscCall(MatGetNonzeroColumnsLocal_Private(PT,&ptap->cis0)); 205 PetscCall(ISLocalToGlobalMappingCreateIS(ptap->cis0,&Ccl2g)); 206 PetscCall(MatDestroy(&PT)); 207 208 Crl2g = NULL; 209 if (rl2g != cl2g) { /* unsymmetric A mapping */ 210 PetscBool same,lsame = PETSC_FALSE; 211 PetscInt N1,ibs1; 212 213 PetscCall(ISLocalToGlobalMappingGetSize(rl2g,&N1)); 214 PetscCall(ISLocalToGlobalMappingGetBlockSize(rl2g,&ibs1)); 215 PetscCall(ISLocalToGlobalMappingGetBlockIndices(rl2g,&garray)); 216 PetscCall(ISCreateBlock(comm,ibs,N/ibs,garray,PETSC_COPY_VALUES,&ptap->ris1)); 217 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(rl2g,&garray)); 218 if (ibs1 == ibs && N1 == N) { /* check if the l2gmaps are the same */ 219 const PetscInt *i1,*i2; 220 221 PetscCall(ISBlockGetIndices(ptap->ris0,&i1)); 222 PetscCall(ISBlockGetIndices(ptap->ris1,&i2)); 223 PetscCall(PetscArraycmp(i1,i2,N,&lsame)); 224 } 225 PetscCall(MPIU_Allreduce(&lsame,&same,1,MPIU_BOOL,MPI_LAND,comm)); 226 if (same) { 227 PetscCall(ISDestroy(&ptap->ris1)); 228 } else { 229 PetscCall(MatCreateSubMatrix(P,ptap->ris1,NULL,MAT_INITIAL_MATRIX,&PT)); 230 PetscCall(MatGetNonzeroColumnsLocal_Private(PT,&ptap->cis1)); 231 PetscCall(ISLocalToGlobalMappingCreateIS(ptap->cis1,&Crl2g)); 232 PetscCall(MatDestroy(&PT)); 233 } 234 } 235 /* Not sure about this 236 if (!Crl2g) { 237 PetscCall(MatGetBlockSize(C,&ibs)); 238 PetscCall(ISLocalToGlobalMappingSetBlockSize(Ccl2g,ibs)); 239 } 240 */ 241 PetscCall(MatSetLocalToGlobalMapping(C,Crl2g ? Crl2g : Ccl2g,Ccl2g)); 242 PetscCall(ISLocalToGlobalMappingDestroy(&Crl2g)); 243 PetscCall(ISLocalToGlobalMappingDestroy(&Ccl2g)); 244 245 C->ops->ptapnumeric = MatPtAPNumeric_IS_XAIJ; 246 PetscFunctionReturn(0); 247 } 248 249 /* ----------------------------------------- */ 250 static PetscErrorCode MatProductSymbolic_PtAP_IS_XAIJ(Mat C) 251 { 252 Mat_Product *product = C->product; 253 Mat A=product->A,P=product->B; 254 PetscReal fill=product->fill; 255 256 PetscFunctionBegin; 257 PetscCall(MatPtAPSymbolic_IS_XAIJ(A,P,fill,C)); 258 C->ops->productnumeric = MatProductNumeric_PtAP; 259 PetscFunctionReturn(0); 260 } 261 262 static PetscErrorCode MatProductSetFromOptions_IS_XAIJ_PtAP(Mat C) 263 { 264 PetscFunctionBegin; 265 C->ops->productsymbolic = MatProductSymbolic_PtAP_IS_XAIJ; 266 PetscFunctionReturn(0); 267 } 268 269 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_IS_XAIJ(Mat C) 270 { 271 Mat_Product *product = C->product; 272 273 PetscFunctionBegin; 274 if (product->type == MATPRODUCT_PtAP) { 275 PetscCall(MatProductSetFromOptions_IS_XAIJ_PtAP(C)); 276 } 277 PetscFunctionReturn(0); 278 } 279 280 /* ----------------------------------------- */ 281 static PetscErrorCode MatISContainerDestroyFields_Private(void *ptr) 282 { 283 MatISLocalFields lf = (MatISLocalFields)ptr; 284 PetscInt i; 285 286 PetscFunctionBegin; 287 for (i=0;i<lf->nr;i++) { 288 PetscCall(ISDestroy(&lf->rf[i])); 289 } 290 for (i=0;i<lf->nc;i++) { 291 PetscCall(ISDestroy(&lf->cf[i])); 292 } 293 PetscCall(PetscFree2(lf->rf,lf->cf)); 294 PetscCall(PetscFree(lf)); 295 PetscFunctionReturn(0); 296 } 297 298 static PetscErrorCode MatConvert_SeqXAIJ_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat) 299 { 300 Mat B,lB; 301 302 PetscFunctionBegin; 303 if (reuse != MAT_REUSE_MATRIX) { 304 ISLocalToGlobalMapping rl2g,cl2g; 305 PetscInt bs; 306 IS is; 307 308 PetscCall(MatGetBlockSize(A,&bs)); 309 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A),A->rmap->n/bs,0,1,&is)); 310 if (bs > 1) { 311 IS is2; 312 PetscInt i,*aux; 313 314 PetscCall(ISGetLocalSize(is,&i)); 315 PetscCall(ISGetIndices(is,(const PetscInt**)&aux)); 316 PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A),bs,i,aux,PETSC_COPY_VALUES,&is2)); 317 PetscCall(ISRestoreIndices(is,(const PetscInt**)&aux)); 318 PetscCall(ISDestroy(&is)); 319 is = is2; 320 } 321 PetscCall(ISSetIdentity(is)); 322 PetscCall(ISLocalToGlobalMappingCreateIS(is,&rl2g)); 323 PetscCall(ISDestroy(&is)); 324 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A),A->cmap->n/bs,0,1,&is)); 325 if (bs > 1) { 326 IS is2; 327 PetscInt i,*aux; 328 329 PetscCall(ISGetLocalSize(is,&i)); 330 PetscCall(ISGetIndices(is,(const PetscInt**)&aux)); 331 PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A),bs,i,aux,PETSC_COPY_VALUES,&is2)); 332 PetscCall(ISRestoreIndices(is,(const PetscInt**)&aux)); 333 PetscCall(ISDestroy(&is)); 334 is = is2; 335 } 336 PetscCall(ISSetIdentity(is)); 337 PetscCall(ISLocalToGlobalMappingCreateIS(is,&cl2g)); 338 PetscCall(ISDestroy(&is)); 339 PetscCall(MatCreateIS(PetscObjectComm((PetscObject)A),bs,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N,rl2g,cl2g,&B)); 340 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 341 PetscCall(ISLocalToGlobalMappingDestroy(&cl2g)); 342 PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&lB)); 343 if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 344 } else { 345 B = *newmat; 346 PetscCall(PetscObjectReference((PetscObject)A)); 347 lB = A; 348 } 349 PetscCall(MatISSetLocalMat(B,lB)); 350 PetscCall(MatDestroy(&lB)); 351 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 352 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 353 if (reuse == MAT_INPLACE_MATRIX) { 354 PetscCall(MatHeaderReplace(A,&B)); 355 } 356 PetscFunctionReturn(0); 357 } 358 359 static PetscErrorCode MatISScaleDisassembling_Private(Mat A) 360 { 361 Mat_IS *matis = (Mat_IS*)(A->data); 362 PetscScalar *aa; 363 const PetscInt *ii,*jj; 364 PetscInt i,n,m; 365 PetscInt *ecount,**eneighs; 366 PetscBool flg; 367 368 PetscFunctionBegin; 369 PetscCall(MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&m,&ii,&jj,&flg)); 370 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get IJ structure"); 371 PetscCall(ISLocalToGlobalMappingGetNodeInfo(matis->rmapping,&n,&ecount,&eneighs)); 372 PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %" PetscInt_FMT " != %" PetscInt_FMT,m,n); 373 PetscCall(MatSeqAIJGetArray(matis->A,&aa)); 374 for (i=0;i<n;i++) { 375 if (ecount[i] > 1) { 376 PetscInt j; 377 378 for (j=ii[i];j<ii[i+1];j++) { 379 PetscInt i2 = jj[j],p,p2; 380 PetscReal scal = 0.0; 381 382 for (p=0;p<ecount[i];p++) { 383 for (p2=0;p2<ecount[i2];p2++) { 384 if (eneighs[i][p] == eneighs[i2][p2]) { scal += 1.0; break; } 385 } 386 } 387 if (scal) aa[j] /= scal; 388 } 389 } 390 } 391 PetscCall(ISLocalToGlobalMappingRestoreNodeInfo(matis->rmapping,&n,&ecount,&eneighs)); 392 PetscCall(MatSeqAIJRestoreArray(matis->A,&aa)); 393 PetscCall(MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&m,&ii,&jj,&flg)); 394 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot restore IJ structure"); 395 PetscFunctionReturn(0); 396 } 397 398 typedef enum {MAT_IS_DISASSEMBLE_L2G_NATURAL,MAT_IS_DISASSEMBLE_L2G_MAT, MAT_IS_DISASSEMBLE_L2G_ND} MatISDisassemblel2gType; 399 400 static PetscErrorCode MatMPIXAIJComputeLocalToGlobalMapping_Private(Mat A, ISLocalToGlobalMapping *l2g) 401 { 402 Mat Ad,Ao; 403 IS is,ndmap,ndsub; 404 MPI_Comm comm; 405 const PetscInt *garray,*ndmapi; 406 PetscInt bs,i,cnt,nl,*ncount,*ndmapc; 407 PetscBool ismpiaij,ismpibaij; 408 const char *const MatISDisassemblel2gTypes[] = {"NATURAL","MAT","ND","MatISDisassemblel2gType","MAT_IS_DISASSEMBLE_L2G_",NULL}; 409 MatISDisassemblel2gType mode = MAT_IS_DISASSEMBLE_L2G_NATURAL; 410 MatPartitioning part; 411 PetscSF sf; 412 413 PetscFunctionBegin; 414 PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatIS l2g disassembling options","Mat"); 415 PetscCall(PetscOptionsEnum("-mat_is_disassemble_l2g_type","Type of local-to-global mapping to be used for disassembling","MatISDisassemblel2gType",MatISDisassemblel2gTypes,(PetscEnum)mode,(PetscEnum*)&mode,NULL)); 416 PetscOptionsEnd(); 417 if (mode == MAT_IS_DISASSEMBLE_L2G_MAT) { 418 PetscCall(MatGetLocalToGlobalMapping(A,l2g,NULL)); 419 PetscFunctionReturn(0); 420 } 421 PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 422 PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATMPIAIJ ,&ismpiaij)); 423 PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATMPIBAIJ,&ismpibaij)); 424 PetscCall(MatGetBlockSize(A,&bs)); 425 switch (mode) { 426 case MAT_IS_DISASSEMBLE_L2G_ND: 427 PetscCall(MatPartitioningCreate(comm,&part)); 428 PetscCall(MatPartitioningSetAdjacency(part,A)); 429 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part,((PetscObject)A)->prefix)); 430 PetscCall(MatPartitioningSetFromOptions(part)); 431 PetscCall(MatPartitioningApplyND(part,&ndmap)); 432 PetscCall(MatPartitioningDestroy(&part)); 433 PetscCall(ISBuildTwoSided(ndmap,NULL,&ndsub)); 434 PetscCall(MatMPIAIJSetUseScalableIncreaseOverlap(A,PETSC_TRUE)); 435 PetscCall(MatIncreaseOverlap(A,1,&ndsub,1)); 436 PetscCall(ISLocalToGlobalMappingCreateIS(ndsub,l2g)); 437 438 /* it may happen that a separator node is not properly shared */ 439 PetscCall(ISLocalToGlobalMappingGetNodeInfo(*l2g,&nl,&ncount,NULL)); 440 PetscCall(PetscSFCreate(comm,&sf)); 441 PetscCall(ISLocalToGlobalMappingGetIndices(*l2g,&garray)); 442 PetscCall(PetscSFSetGraphLayout(sf,A->rmap,nl,NULL,PETSC_OWN_POINTER,garray)); 443 PetscCall(ISLocalToGlobalMappingRestoreIndices(*l2g,&garray)); 444 PetscCall(PetscCalloc1(A->rmap->n,&ndmapc)); 445 PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ncount,ndmapc,MPI_REPLACE)); 446 PetscCall(PetscSFReduceEnd(sf,MPIU_INT,ncount,ndmapc,MPI_REPLACE)); 447 PetscCall(ISLocalToGlobalMappingRestoreNodeInfo(*l2g,NULL,&ncount,NULL)); 448 PetscCall(ISGetIndices(ndmap,&ndmapi)); 449 for (i = 0, cnt = 0; i < A->rmap->n; i++) 450 if (ndmapi[i] < 0 && ndmapc[i] < 2) 451 cnt++; 452 453 PetscCall(MPIU_Allreduce(&cnt,&i,1,MPIU_INT,MPI_MAX,comm)); 454 if (i) { /* we detected isolated separator nodes */ 455 Mat A2,A3; 456 IS *workis,is2; 457 PetscScalar *vals; 458 PetscInt gcnt = i,*dnz,*onz,j,*lndmapi; 459 ISLocalToGlobalMapping ll2g; 460 PetscBool flg; 461 const PetscInt *ii,*jj; 462 463 /* communicate global id of separators */ 464 MatPreallocateBegin(comm,A->rmap->n,A->cmap->n,dnz,onz); 465 for (i = 0, cnt = 0; i < A->rmap->n; i++) 466 dnz[i] = ndmapi[i] < 0 ? i + A->rmap->rstart : -1; 467 468 PetscCall(PetscMalloc1(nl,&lndmapi)); 469 PetscCall(PetscSFBcastBegin(sf,MPIU_INT,dnz,lndmapi,MPI_REPLACE)); 470 471 /* compute adjacency of isolated separators node */ 472 PetscCall(PetscMalloc1(gcnt,&workis)); 473 for (i = 0, cnt = 0; i < A->rmap->n; i++) { 474 if (ndmapi[i] < 0 && ndmapc[i] < 2) { 475 PetscCall(ISCreateStride(comm,1,i+A->rmap->rstart,1,&workis[cnt++])); 476 } 477 } 478 for (i = cnt; i < gcnt; i++) { 479 PetscCall(ISCreateStride(comm,0,0,1,&workis[i])); 480 } 481 for (i = 0; i < gcnt; i++) { 482 PetscCall(PetscObjectSetName((PetscObject)workis[i],"ISOLATED")); 483 PetscCall(ISViewFromOptions(workis[i],NULL,"-view_isolated_separators")); 484 } 485 486 /* no communications since all the ISes correspond to locally owned rows */ 487 PetscCall(MatIncreaseOverlap(A,gcnt,workis,1)); 488 489 /* end communicate global id of separators */ 490 PetscCall(PetscSFBcastEnd(sf,MPIU_INT,dnz,lndmapi,MPI_REPLACE)); 491 492 /* communicate new layers : create a matrix and transpose it */ 493 PetscCall(PetscArrayzero(dnz,A->rmap->n)); 494 PetscCall(PetscArrayzero(onz,A->rmap->n)); 495 for (i = 0, j = 0; i < A->rmap->n; i++) { 496 if (ndmapi[i] < 0 && ndmapc[i] < 2) { 497 const PetscInt* idxs; 498 PetscInt s; 499 500 PetscCall(ISGetLocalSize(workis[j],&s)); 501 PetscCall(ISGetIndices(workis[j],&idxs)); 502 PetscCall(MatPreallocateSet(i+A->rmap->rstart,s,idxs,dnz,onz)); 503 j++; 504 } 505 } 506 PetscCheck(j == cnt,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected local count %" PetscInt_FMT " != %" PetscInt_FMT,j,cnt); 507 508 for (i = 0; i < gcnt; i++) { 509 PetscCall(PetscObjectSetName((PetscObject)workis[i],"EXTENDED")); 510 PetscCall(ISViewFromOptions(workis[i],NULL,"-view_isolated_separators")); 511 } 512 513 for (i = 0, j = 0; i < A->rmap->n; i++) j = PetscMax(j,dnz[i]+onz[i]); 514 PetscCall(PetscMalloc1(j,&vals)); 515 for (i = 0; i < j; i++) vals[i] = 1.0; 516 517 PetscCall(MatCreate(comm,&A2)); 518 PetscCall(MatSetType(A2,MATMPIAIJ)); 519 PetscCall(MatSetSizes(A2,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 520 PetscCall(MatMPIAIJSetPreallocation(A2,0,dnz,0,onz)); 521 PetscCall(MatSetOption(A2,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 522 for (i = 0, j = 0; i < A2->rmap->n; i++) { 523 PetscInt row = i+A2->rmap->rstart,s = dnz[i] + onz[i]; 524 const PetscInt* idxs; 525 526 if (s) { 527 PetscCall(ISGetIndices(workis[j],&idxs)); 528 PetscCall(MatSetValues(A2,1,&row,s,idxs,vals,INSERT_VALUES)); 529 PetscCall(ISRestoreIndices(workis[j],&idxs)); 530 j++; 531 } 532 } 533 PetscCheck(j == cnt,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected local count %" PetscInt_FMT " != %" PetscInt_FMT,j,cnt); 534 PetscCall(PetscFree(vals)); 535 PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 536 PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 537 PetscCall(MatTranspose(A2,MAT_INPLACE_MATRIX,&A2)); 538 539 /* extract submatrix corresponding to the coupling "owned separators" x "isolated separators" */ 540 for (i = 0, j = 0; i < nl; i++) 541 if (lndmapi[i] >= 0) lndmapi[j++] = lndmapi[i]; 542 PetscCall(ISCreateGeneral(comm,j,lndmapi,PETSC_USE_POINTER,&is)); 543 PetscCall(MatMPIAIJGetLocalMatCondensed(A2,MAT_INITIAL_MATRIX,&is,NULL,&A3)); 544 PetscCall(ISDestroy(&is)); 545 PetscCall(MatDestroy(&A2)); 546 547 /* extend local to global map to include connected isolated separators */ 548 PetscCall(PetscObjectQuery((PetscObject)A3,"_petsc_GetLocalMatCondensed_iscol",(PetscObject*)&is)); 549 PetscCheck(is,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing column map"); 550 PetscCall(ISLocalToGlobalMappingCreateIS(is,&ll2g)); 551 PetscCall(MatGetRowIJ(A3,0,PETSC_FALSE,PETSC_FALSE,&i,&ii,&jj,&flg)); 552 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get IJ structure"); 553 PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ii[i],jj,PETSC_COPY_VALUES,&is)); 554 PetscCall(MatRestoreRowIJ(A3,0,PETSC_FALSE,PETSC_FALSE,&i,&ii,&jj,&flg)); 555 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get IJ structure"); 556 PetscCall(ISLocalToGlobalMappingApplyIS(ll2g,is,&is2)); 557 PetscCall(ISDestroy(&is)); 558 PetscCall(ISLocalToGlobalMappingDestroy(&ll2g)); 559 560 /* add new nodes to the local-to-global map */ 561 PetscCall(ISLocalToGlobalMappingDestroy(l2g)); 562 PetscCall(ISExpand(ndsub,is2,&is)); 563 PetscCall(ISDestroy(&is2)); 564 PetscCall(ISLocalToGlobalMappingCreateIS(is,l2g)); 565 PetscCall(ISDestroy(&is)); 566 567 PetscCall(MatDestroy(&A3)); 568 PetscCall(PetscFree(lndmapi)); 569 MatPreallocateEnd(dnz,onz); 570 for (i = 0; i < gcnt; i++) { 571 PetscCall(ISDestroy(&workis[i])); 572 } 573 PetscCall(PetscFree(workis)); 574 } 575 PetscCall(ISRestoreIndices(ndmap,&ndmapi)); 576 PetscCall(PetscSFDestroy(&sf)); 577 PetscCall(PetscFree(ndmapc)); 578 PetscCall(ISDestroy(&ndmap)); 579 PetscCall(ISDestroy(&ndsub)); 580 PetscCall(ISLocalToGlobalMappingSetBlockSize(*l2g,bs)); 581 PetscCall(ISLocalToGlobalMappingViewFromOptions(*l2g,NULL,"-matis_nd_l2g_view")); 582 break; 583 case MAT_IS_DISASSEMBLE_L2G_NATURAL: 584 if (ismpiaij) { 585 PetscCall(MatMPIAIJGetSeqAIJ(A,&Ad,&Ao,&garray)); 586 } else if (ismpibaij) { 587 PetscCall(MatMPIBAIJGetSeqBAIJ(A,&Ad,&Ao,&garray)); 588 } else SETERRQ(comm,PETSC_ERR_SUP,"Type %s",((PetscObject)A)->type_name); 589 PetscCheck(garray,comm,PETSC_ERR_ARG_WRONGSTATE,"garray not present"); 590 if (A->rmap->n) { 591 PetscInt dc,oc,stc,*aux; 592 593 PetscCall(MatGetLocalSize(Ad,NULL,&dc)); 594 PetscCall(MatGetLocalSize(Ao,NULL,&oc)); 595 PetscCall(MatGetOwnershipRangeColumn(A,&stc,NULL)); 596 PetscCall(PetscMalloc1((dc+oc)/bs,&aux)); 597 for (i=0; i<dc/bs; i++) aux[i] = i+stc/bs; 598 for (i=0; i<oc/bs; i++) aux[i+dc/bs] = (ismpiaij ? garray[i*bs]/bs : garray[i]); 599 PetscCall(ISCreateBlock(comm,bs,(dc+oc)/bs,aux,PETSC_OWN_POINTER,&is)); 600 } else { 601 PetscCall(ISCreateBlock(comm,1,0,NULL,PETSC_OWN_POINTER,&is)); 602 } 603 PetscCall(ISLocalToGlobalMappingCreateIS(is,l2g)); 604 PetscCall(ISDestroy(&is)); 605 break; 606 default: 607 SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Unsupported l2g disassembling type %d",mode); 608 } 609 PetscFunctionReturn(0); 610 } 611 612 PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat) 613 { 614 Mat lA,Ad,Ao,B = NULL; 615 ISLocalToGlobalMapping rl2g,cl2g; 616 IS is; 617 MPI_Comm comm; 618 void *ptrs[2]; 619 const char *names[2] = {"_convert_csr_aux","_convert_csr_data"}; 620 const PetscInt *garray; 621 PetscScalar *dd,*od,*aa,*data; 622 const PetscInt *di,*dj,*oi,*oj; 623 const PetscInt *odi,*odj,*ooi,*ooj; 624 PetscInt *aux,*ii,*jj; 625 PetscInt bs,lc,dr,dc,oc,str,stc,nnz,i,jd,jo,cum; 626 PetscBool flg,ismpiaij,ismpibaij,was_inplace = PETSC_FALSE; 627 PetscMPIInt size; 628 629 PetscFunctionBegin; 630 PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 631 PetscCallMPI(MPI_Comm_size(comm,&size)); 632 if (size == 1) { 633 PetscCall(MatConvert_SeqXAIJ_IS(A,type,reuse,newmat)); 634 PetscFunctionReturn(0); 635 } 636 if (reuse != MAT_REUSE_MATRIX && A->cmap->N == A->rmap->N) { 637 PetscCall(MatMPIXAIJComputeLocalToGlobalMapping_Private(A,&rl2g)); 638 PetscCall(MatCreate(comm,&B)); 639 PetscCall(MatSetType(B,MATIS)); 640 PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 641 PetscCall(MatSetLocalToGlobalMapping(B,rl2g,rl2g)); 642 PetscCall(MatGetBlockSize(A,&bs)); 643 PetscCall(MatSetBlockSize(B,bs)); 644 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 645 if (reuse == MAT_INPLACE_MATRIX) was_inplace = PETSC_TRUE; 646 reuse = MAT_REUSE_MATRIX; 647 } 648 if (reuse == MAT_REUSE_MATRIX) { 649 Mat *newlA, lA; 650 IS rows, cols; 651 const PetscInt *ridx, *cidx; 652 PetscInt rbs, cbs, nr, nc; 653 654 if (!B) B = *newmat; 655 PetscCall(MatISGetLocalToGlobalMapping(B,&rl2g,&cl2g)); 656 PetscCall(ISLocalToGlobalMappingGetBlockIndices(rl2g,&ridx)); 657 PetscCall(ISLocalToGlobalMappingGetBlockIndices(cl2g,&cidx)); 658 PetscCall(ISLocalToGlobalMappingGetSize(rl2g,&nr)); 659 PetscCall(ISLocalToGlobalMappingGetSize(cl2g,&nc)); 660 PetscCall(ISLocalToGlobalMappingGetBlockSize(rl2g,&rbs)); 661 PetscCall(ISLocalToGlobalMappingGetBlockSize(cl2g,&cbs)); 662 PetscCall(ISCreateBlock(comm,rbs,nr/rbs,ridx,PETSC_USE_POINTER,&rows)); 663 if (rl2g != cl2g) { 664 PetscCall(ISCreateBlock(comm,cbs,nc/cbs,cidx,PETSC_USE_POINTER,&cols)); 665 } else { 666 PetscCall(PetscObjectReference((PetscObject)rows)); 667 cols = rows; 668 } 669 PetscCall(MatISGetLocalMat(B,&lA)); 670 PetscCall(MatCreateSubMatrices(A,1,&rows,&cols,MAT_INITIAL_MATRIX,&newlA)); 671 PetscCall(MatConvert(newlA[0],MATSEQAIJ,MAT_INPLACE_MATRIX,&newlA[0])); 672 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(rl2g,&ridx)); 673 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(cl2g,&cidx)); 674 PetscCall(ISDestroy(&rows)); 675 PetscCall(ISDestroy(&cols)); 676 if (!lA->preallocated) { /* first time */ 677 PetscCall(MatDuplicate(newlA[0],MAT_COPY_VALUES,&lA)); 678 PetscCall(MatISSetLocalMat(B,lA)); 679 PetscCall(PetscObjectDereference((PetscObject)lA)); 680 } 681 PetscCall(MatCopy(newlA[0],lA,SAME_NONZERO_PATTERN)); 682 PetscCall(MatDestroySubMatrices(1,&newlA)); 683 PetscCall(MatISScaleDisassembling_Private(B)); 684 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 685 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 686 if (was_inplace) PetscCall(MatHeaderReplace(A,&B)); 687 else *newmat = B; 688 PetscFunctionReturn(0); 689 } 690 /* rectangular case, just compress out the column space */ 691 PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATMPIAIJ ,&ismpiaij)); 692 PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATMPIBAIJ,&ismpibaij)); 693 if (ismpiaij) { 694 bs = 1; 695 PetscCall(MatMPIAIJGetSeqAIJ(A,&Ad,&Ao,&garray)); 696 } else if (ismpibaij) { 697 PetscCall(MatGetBlockSize(A,&bs)); 698 PetscCall(MatMPIBAIJGetSeqBAIJ(A,&Ad,&Ao,&garray)); 699 PetscCall(MatConvert(Ad,MATSEQAIJ,MAT_INITIAL_MATRIX,&Ad)); 700 PetscCall(MatConvert(Ao,MATSEQAIJ,MAT_INITIAL_MATRIX,&Ao)); 701 } else SETERRQ(comm,PETSC_ERR_SUP,"Type %s",((PetscObject)A)->type_name); 702 PetscCall(MatSeqAIJGetArray(Ad,&dd)); 703 PetscCall(MatSeqAIJGetArray(Ao,&od)); 704 PetscCheck(garray,comm,PETSC_ERR_ARG_WRONGSTATE,"garray not present"); 705 706 /* access relevant information from MPIAIJ */ 707 PetscCall(MatGetOwnershipRange(A,&str,NULL)); 708 PetscCall(MatGetOwnershipRangeColumn(A,&stc,NULL)); 709 PetscCall(MatGetLocalSize(A,&dr,&dc)); 710 PetscCall(MatGetLocalSize(Ao,NULL,&oc)); 711 PetscCall(MatGetRowIJ(Ad,0,PETSC_FALSE,PETSC_FALSE,&i,&di,&dj,&flg)); 712 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get IJ structure"); 713 PetscCall(MatGetRowIJ(Ao,0,PETSC_FALSE,PETSC_FALSE,&i,&oi,&oj,&flg)); 714 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get IJ structure"); 715 nnz = di[dr] + oi[dr]; 716 /* store original pointers to be restored later */ 717 odi = di; odj = dj; ooi = oi; ooj = oj; 718 719 /* generate l2g maps for rows and cols */ 720 PetscCall(ISCreateStride(comm,dr/bs,str/bs,1,&is)); 721 if (bs > 1) { 722 IS is2; 723 724 PetscCall(ISGetLocalSize(is,&i)); 725 PetscCall(ISGetIndices(is,(const PetscInt**)&aux)); 726 PetscCall(ISCreateBlock(comm,bs,i,aux,PETSC_COPY_VALUES,&is2)); 727 PetscCall(ISRestoreIndices(is,(const PetscInt**)&aux)); 728 PetscCall(ISDestroy(&is)); 729 is = is2; 730 } 731 PetscCall(ISLocalToGlobalMappingCreateIS(is,&rl2g)); 732 PetscCall(ISDestroy(&is)); 733 if (dr) { 734 PetscCall(PetscMalloc1((dc+oc)/bs,&aux)); 735 for (i=0; i<dc/bs; i++) aux[i] = i+stc/bs; 736 for (i=0; i<oc/bs; i++) aux[i+dc/bs] = garray[i]; 737 PetscCall(ISCreateBlock(comm,bs,(dc+oc)/bs,aux,PETSC_OWN_POINTER,&is)); 738 lc = dc+oc; 739 } else { 740 PetscCall(ISCreateBlock(comm,bs,0,NULL,PETSC_OWN_POINTER,&is)); 741 lc = 0; 742 } 743 PetscCall(ISLocalToGlobalMappingCreateIS(is,&cl2g)); 744 PetscCall(ISDestroy(&is)); 745 746 /* create MATIS object */ 747 PetscCall(MatCreate(comm,&B)); 748 PetscCall(MatSetSizes(B,dr,dc,PETSC_DECIDE,PETSC_DECIDE)); 749 PetscCall(MatSetType(B,MATIS)); 750 PetscCall(MatSetBlockSize(B,bs)); 751 PetscCall(MatSetLocalToGlobalMapping(B,rl2g,cl2g)); 752 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 753 PetscCall(ISLocalToGlobalMappingDestroy(&cl2g)); 754 755 /* merge local matrices */ 756 PetscCall(PetscMalloc1(nnz+dr+1,&aux)); 757 PetscCall(PetscMalloc1(nnz,&data)); 758 ii = aux; 759 jj = aux+dr+1; 760 aa = data; 761 *ii = *(di++) + *(oi++); 762 for (jd=0,jo=0,cum=0;*ii<nnz;cum++) 763 { 764 for (;jd<*di;jd++) { *jj++ = *dj++; *aa++ = *dd++; } 765 for (;jo<*oi;jo++) { *jj++ = *oj++ + dc; *aa++ = *od++; } 766 *(++ii) = *(di++) + *(oi++); 767 } 768 for (;cum<dr;cum++) *(++ii) = nnz; 769 770 PetscCall(MatRestoreRowIJ(Ad,0,PETSC_FALSE,PETSC_FALSE,&i,&odi,&odj,&flg)); 771 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot restore IJ structure"); 772 PetscCall(MatRestoreRowIJ(Ao,0,PETSC_FALSE,PETSC_FALSE,&i,&ooi,&ooj,&flg)); 773 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot restore IJ structure"); 774 PetscCall(MatSeqAIJRestoreArray(Ad,&dd)); 775 PetscCall(MatSeqAIJRestoreArray(Ao,&od)); 776 777 ii = aux; 778 jj = aux+dr+1; 779 aa = data; 780 PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,lc,ii,jj,aa,&lA)); 781 782 /* create containers to destroy the data */ 783 ptrs[0] = aux; 784 ptrs[1] = data; 785 for (i=0; i<2; i++) { 786 PetscContainer c; 787 788 PetscCall(PetscContainerCreate(PETSC_COMM_SELF,&c)); 789 PetscCall(PetscContainerSetPointer(c,ptrs[i])); 790 PetscCall(PetscContainerSetUserDestroy(c,PetscContainerUserDestroyDefault)); 791 PetscCall(PetscObjectCompose((PetscObject)lA,names[i],(PetscObject)c)); 792 PetscCall(PetscContainerDestroy(&c)); 793 } 794 if (ismpibaij) { /* destroy converted local matrices */ 795 PetscCall(MatDestroy(&Ad)); 796 PetscCall(MatDestroy(&Ao)); 797 } 798 799 /* finalize matrix */ 800 PetscCall(MatISSetLocalMat(B,lA)); 801 PetscCall(MatDestroy(&lA)); 802 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 803 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 804 if (reuse == MAT_INPLACE_MATRIX) { 805 PetscCall(MatHeaderReplace(A,&B)); 806 } else *newmat = B; 807 PetscFunctionReturn(0); 808 } 809 810 PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat) 811 { 812 Mat **nest,*snest,**rnest,lA,B; 813 IS *iscol,*isrow,*islrow,*islcol; 814 ISLocalToGlobalMapping rl2g,cl2g; 815 MPI_Comm comm; 816 PetscInt *lr,*lc,*l2gidxs; 817 PetscInt i,j,nr,nc,rbs,cbs; 818 PetscBool convert,lreuse,*istrans; 819 820 PetscFunctionBegin; 821 PetscCall(MatNestGetSubMats(A,&nr,&nc,&nest)); 822 lreuse = PETSC_FALSE; 823 rnest = NULL; 824 if (reuse == MAT_REUSE_MATRIX) { 825 PetscBool ismatis,isnest; 826 827 PetscCall(PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis)); 828 PetscCheck(ismatis,PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type_name); 829 PetscCall(MatISGetLocalMat(*newmat,&lA)); 830 PetscCall(PetscObjectTypeCompare((PetscObject)lA,MATNEST,&isnest)); 831 if (isnest) { 832 PetscCall(MatNestGetSubMats(lA,&i,&j,&rnest)); 833 lreuse = (PetscBool)(i == nr && j == nc); 834 if (!lreuse) rnest = NULL; 835 } 836 } 837 PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 838 PetscCall(PetscCalloc2(nr,&lr,nc,&lc)); 839 PetscCall(PetscCalloc6(nr,&isrow,nc,&iscol,nr,&islrow,nc,&islcol,nr*nc,&snest,nr*nc,&istrans)); 840 PetscCall(MatNestGetISs(A,isrow,iscol)); 841 for (i=0;i<nr;i++) { 842 for (j=0;j<nc;j++) { 843 PetscBool ismatis; 844 PetscInt l1,l2,lb1,lb2,ij=i*nc+j; 845 846 /* Null matrix pointers are allowed in MATNEST */ 847 if (!nest[i][j]) continue; 848 849 /* Nested matrices should be of type MATIS */ 850 PetscCall(PetscObjectTypeCompare((PetscObject)nest[i][j],MATTRANSPOSEMAT,&istrans[ij])); 851 if (istrans[ij]) { 852 Mat T,lT; 853 PetscCall(MatTransposeGetMat(nest[i][j],&T)); 854 PetscCall(PetscObjectTypeCompare((PetscObject)T,MATIS,&ismatis)); 855 PetscCheck(ismatis,comm,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") (transposed) is not of type MATIS",i,j); 856 PetscCall(MatISGetLocalMat(T,&lT)); 857 PetscCall(MatCreateTranspose(lT,&snest[ij])); 858 } else { 859 PetscCall(PetscObjectTypeCompare((PetscObject)nest[i][j],MATIS,&ismatis)); 860 PetscCheck(ismatis,comm,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") is not of type MATIS",i,j); 861 PetscCall(MatISGetLocalMat(nest[i][j],&snest[ij])); 862 } 863 864 /* Check compatibility of local sizes */ 865 PetscCall(MatGetSize(snest[ij],&l1,&l2)); 866 PetscCall(MatGetBlockSizes(snest[ij],&lb1,&lb2)); 867 if (!l1 || !l2) continue; 868 PetscCheckFalse(lr[i] && l1 != lr[i],PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid local size %" PetscInt_FMT " != %" PetscInt_FMT,i,j,lr[i],l1); 869 PetscCheckFalse(lc[j] && l2 != lc[j],PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid local size %" PetscInt_FMT " != %" PetscInt_FMT,i,j,lc[j],l2); 870 lr[i] = l1; 871 lc[j] = l2; 872 873 /* check compatibilty for local matrix reusage */ 874 if (rnest && !rnest[i][j] != !snest[ij]) lreuse = PETSC_FALSE; 875 } 876 } 877 878 if (PetscDefined (USE_DEBUG)) { 879 /* Check compatibility of l2g maps for rows */ 880 for (i=0;i<nr;i++) { 881 rl2g = NULL; 882 for (j=0;j<nc;j++) { 883 PetscInt n1,n2; 884 885 if (!nest[i][j]) continue; 886 if (istrans[i*nc+j]) { 887 Mat T; 888 889 PetscCall(MatTransposeGetMat(nest[i][j],&T)); 890 PetscCall(MatISGetLocalToGlobalMapping(T,NULL,&cl2g)); 891 } else { 892 PetscCall(MatISGetLocalToGlobalMapping(nest[i][j],&cl2g,NULL)); 893 } 894 PetscCall(ISLocalToGlobalMappingGetSize(cl2g,&n1)); 895 if (!n1) continue; 896 if (!rl2g) { 897 rl2g = cl2g; 898 } else { 899 const PetscInt *idxs1,*idxs2; 900 PetscBool same; 901 902 PetscCall(ISLocalToGlobalMappingGetSize(rl2g,&n2)); 903 PetscCheck(n1 == n2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid row l2gmap size %" PetscInt_FMT " != %" PetscInt_FMT,i,j,n1,n2); 904 PetscCall(ISLocalToGlobalMappingGetIndices(cl2g,&idxs1)); 905 PetscCall(ISLocalToGlobalMappingGetIndices(rl2g,&idxs2)); 906 PetscCall(PetscArraycmp(idxs1,idxs2,n1,&same)); 907 PetscCall(ISLocalToGlobalMappingRestoreIndices(cl2g,&idxs1)); 908 PetscCall(ISLocalToGlobalMappingRestoreIndices(rl2g,&idxs2)); 909 PetscCheck(same,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid row l2gmap",i,j); 910 } 911 } 912 } 913 /* Check compatibility of l2g maps for columns */ 914 for (i=0;i<nc;i++) { 915 rl2g = NULL; 916 for (j=0;j<nr;j++) { 917 PetscInt n1,n2; 918 919 if (!nest[j][i]) continue; 920 if (istrans[j*nc+i]) { 921 Mat T; 922 923 PetscCall(MatTransposeGetMat(nest[j][i],&T)); 924 PetscCall(MatISGetLocalToGlobalMapping(T,&cl2g,NULL)); 925 } else { 926 PetscCall(MatISGetLocalToGlobalMapping(nest[j][i],NULL,&cl2g)); 927 } 928 PetscCall(ISLocalToGlobalMappingGetSize(cl2g,&n1)); 929 if (!n1) continue; 930 if (!rl2g) { 931 rl2g = cl2g; 932 } else { 933 const PetscInt *idxs1,*idxs2; 934 PetscBool same; 935 936 PetscCall(ISLocalToGlobalMappingGetSize(rl2g,&n2)); 937 PetscCheck(n1 == n2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid column l2gmap size %" PetscInt_FMT " != %" PetscInt_FMT,j,i,n1,n2); 938 PetscCall(ISLocalToGlobalMappingGetIndices(cl2g,&idxs1)); 939 PetscCall(ISLocalToGlobalMappingGetIndices(rl2g,&idxs2)); 940 PetscCall(PetscArraycmp(idxs1,idxs2,n1,&same)); 941 PetscCall(ISLocalToGlobalMappingRestoreIndices(cl2g,&idxs1)); 942 PetscCall(ISLocalToGlobalMappingRestoreIndices(rl2g,&idxs2)); 943 PetscCheck(same,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid column l2gmap",j,i); 944 } 945 } 946 } 947 } 948 949 B = NULL; 950 if (reuse != MAT_REUSE_MATRIX) { 951 PetscInt stl; 952 953 /* Create l2g map for the rows of the new matrix and index sets for the local MATNEST */ 954 for (i=0,stl=0;i<nr;i++) stl += lr[i]; 955 PetscCall(PetscMalloc1(stl,&l2gidxs)); 956 for (i=0,stl=0;i<nr;i++) { 957 Mat usedmat; 958 Mat_IS *matis; 959 const PetscInt *idxs; 960 961 /* local IS for local NEST */ 962 PetscCall(ISCreateStride(PETSC_COMM_SELF,lr[i],stl,1,&islrow[i])); 963 964 /* l2gmap */ 965 j = 0; 966 usedmat = nest[i][j]; 967 while (!usedmat && j < nc-1) usedmat = nest[i][++j]; 968 PetscCheck(usedmat,comm,PETSC_ERR_SUP,"Cannot find valid row mat"); 969 970 if (istrans[i*nc+j]) { 971 Mat T; 972 PetscCall(MatTransposeGetMat(usedmat,&T)); 973 usedmat = T; 974 } 975 matis = (Mat_IS*)(usedmat->data); 976 PetscCall(ISGetIndices(isrow[i],&idxs)); 977 if (istrans[i*nc+j]) { 978 PetscCall(PetscSFBcastBegin(matis->csf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 979 PetscCall(PetscSFBcastEnd(matis->csf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 980 } else { 981 PetscCall(PetscSFBcastBegin(matis->sf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 982 PetscCall(PetscSFBcastEnd(matis->sf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 983 } 984 PetscCall(ISRestoreIndices(isrow[i],&idxs)); 985 stl += lr[i]; 986 } 987 PetscCall(ISLocalToGlobalMappingCreate(comm,1,stl,l2gidxs,PETSC_OWN_POINTER,&rl2g)); 988 989 /* Create l2g map for columns of the new matrix and index sets for the local MATNEST */ 990 for (i=0,stl=0;i<nc;i++) stl += lc[i]; 991 PetscCall(PetscMalloc1(stl,&l2gidxs)); 992 for (i=0,stl=0;i<nc;i++) { 993 Mat usedmat; 994 Mat_IS *matis; 995 const PetscInt *idxs; 996 997 /* local IS for local NEST */ 998 PetscCall(ISCreateStride(PETSC_COMM_SELF,lc[i],stl,1,&islcol[i])); 999 1000 /* l2gmap */ 1001 j = 0; 1002 usedmat = nest[j][i]; 1003 while (!usedmat && j < nr-1) usedmat = nest[++j][i]; 1004 PetscCheck(usedmat,comm,PETSC_ERR_SUP,"Cannot find valid column mat"); 1005 if (istrans[j*nc+i]) { 1006 Mat T; 1007 PetscCall(MatTransposeGetMat(usedmat,&T)); 1008 usedmat = T; 1009 } 1010 matis = (Mat_IS*)(usedmat->data); 1011 PetscCall(ISGetIndices(iscol[i],&idxs)); 1012 if (istrans[j*nc+i]) { 1013 PetscCall(PetscSFBcastBegin(matis->sf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 1014 PetscCall(PetscSFBcastEnd(matis->sf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 1015 } else { 1016 PetscCall(PetscSFBcastBegin(matis->csf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 1017 PetscCall(PetscSFBcastEnd(matis->csf,MPIU_INT,idxs,l2gidxs+stl,MPI_REPLACE)); 1018 } 1019 PetscCall(ISRestoreIndices(iscol[i],&idxs)); 1020 stl += lc[i]; 1021 } 1022 PetscCall(ISLocalToGlobalMappingCreate(comm,1,stl,l2gidxs,PETSC_OWN_POINTER,&cl2g)); 1023 1024 /* Create MATIS */ 1025 PetscCall(MatCreate(comm,&B)); 1026 PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 1027 PetscCall(MatGetBlockSizes(A,&rbs,&cbs)); 1028 PetscCall(MatSetBlockSizes(B,rbs,cbs)); 1029 PetscCall(MatSetType(B,MATIS)); 1030 PetscCall(MatISSetLocalMatType(B,MATNEST)); 1031 { /* hack : avoid setup of scatters */ 1032 Mat_IS *matis = (Mat_IS*)(B->data); 1033 matis->islocalref = PETSC_TRUE; 1034 } 1035 PetscCall(MatSetLocalToGlobalMapping(B,rl2g,cl2g)); 1036 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 1037 PetscCall(ISLocalToGlobalMappingDestroy(&cl2g)); 1038 PetscCall(MatCreateNest(PETSC_COMM_SELF,nr,islrow,nc,islcol,snest,&lA)); 1039 PetscCall(MatNestSetVecType(lA,VECNEST)); 1040 for (i=0;i<nr*nc;i++) { 1041 if (istrans[i]) { 1042 PetscCall(MatDestroy(&snest[i])); 1043 } 1044 } 1045 PetscCall(MatISSetLocalMat(B,lA)); 1046 PetscCall(MatDestroy(&lA)); 1047 { /* hack : setup of scatters done here */ 1048 Mat_IS *matis = (Mat_IS*)(B->data); 1049 1050 matis->islocalref = PETSC_FALSE; 1051 PetscCall(MatISSetUpScatters_Private(B)); 1052 } 1053 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1054 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1055 if (reuse == MAT_INPLACE_MATRIX) { 1056 PetscCall(MatHeaderReplace(A,&B)); 1057 } else { 1058 *newmat = B; 1059 } 1060 } else { 1061 if (lreuse) { 1062 PetscCall(MatISGetLocalMat(*newmat,&lA)); 1063 for (i=0;i<nr;i++) { 1064 for (j=0;j<nc;j++) { 1065 if (snest[i*nc+j]) { 1066 PetscCall(MatNestSetSubMat(lA,i,j,snest[i*nc+j])); 1067 if (istrans[i*nc+j]) { 1068 PetscCall(MatDestroy(&snest[i*nc+j])); 1069 } 1070 } 1071 } 1072 } 1073 } else { 1074 PetscInt stl; 1075 for (i=0,stl=0;i<nr;i++) { 1076 PetscCall(ISCreateStride(PETSC_COMM_SELF,lr[i],stl,1,&islrow[i])); 1077 stl += lr[i]; 1078 } 1079 for (i=0,stl=0;i<nc;i++) { 1080 PetscCall(ISCreateStride(PETSC_COMM_SELF,lc[i],stl,1,&islcol[i])); 1081 stl += lc[i]; 1082 } 1083 PetscCall(MatCreateNest(PETSC_COMM_SELF,nr,islrow,nc,islcol,snest,&lA)); 1084 for (i=0;i<nr*nc;i++) { 1085 if (istrans[i]) { 1086 PetscCall(MatDestroy(&snest[i])); 1087 } 1088 } 1089 PetscCall(MatISSetLocalMat(*newmat,lA)); 1090 PetscCall(MatDestroy(&lA)); 1091 } 1092 PetscCall(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY)); 1093 PetscCall(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY)); 1094 } 1095 1096 /* Create local matrix in MATNEST format */ 1097 convert = PETSC_FALSE; 1098 PetscCall(PetscOptionsGetBool(NULL,((PetscObject)A)->prefix,"-matis_convert_local_nest",&convert,NULL)); 1099 if (convert) { 1100 Mat M; 1101 MatISLocalFields lf; 1102 PetscContainer c; 1103 1104 PetscCall(MatISGetLocalMat(*newmat,&lA)); 1105 PetscCall(MatConvert(lA,MATAIJ,MAT_INITIAL_MATRIX,&M)); 1106 PetscCall(MatISSetLocalMat(*newmat,M)); 1107 PetscCall(MatDestroy(&M)); 1108 1109 /* attach local fields to the matrix */ 1110 PetscCall(PetscNew(&lf)); 1111 PetscCall(PetscMalloc2(nr,&lf->rf,nc,&lf->cf)); 1112 for (i=0;i<nr;i++) { 1113 PetscInt n,st; 1114 1115 PetscCall(ISGetLocalSize(islrow[i],&n)); 1116 PetscCall(ISStrideGetInfo(islrow[i],&st,NULL)); 1117 PetscCall(ISCreateStride(comm,n,st,1,&lf->rf[i])); 1118 } 1119 for (i=0;i<nc;i++) { 1120 PetscInt n,st; 1121 1122 PetscCall(ISGetLocalSize(islcol[i],&n)); 1123 PetscCall(ISStrideGetInfo(islcol[i],&st,NULL)); 1124 PetscCall(ISCreateStride(comm,n,st,1,&lf->cf[i])); 1125 } 1126 lf->nr = nr; 1127 lf->nc = nc; 1128 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)(*newmat)),&c)); 1129 PetscCall(PetscContainerSetPointer(c,lf)); 1130 PetscCall(PetscContainerSetUserDestroy(c,MatISContainerDestroyFields_Private)); 1131 PetscCall(PetscObjectCompose((PetscObject)(*newmat),"_convert_nest_lfields",(PetscObject)c)); 1132 PetscCall(PetscContainerDestroy(&c)); 1133 } 1134 1135 /* Free workspace */ 1136 for (i=0;i<nr;i++) { 1137 PetscCall(ISDestroy(&islrow[i])); 1138 } 1139 for (i=0;i<nc;i++) { 1140 PetscCall(ISDestroy(&islcol[i])); 1141 } 1142 PetscCall(PetscFree6(isrow,iscol,islrow,islcol,snest,istrans)); 1143 PetscCall(PetscFree2(lr,lc)); 1144 PetscFunctionReturn(0); 1145 } 1146 1147 static PetscErrorCode MatDiagonalScale_IS(Mat A, Vec l, Vec r) 1148 { 1149 Mat_IS *matis = (Mat_IS*)A->data; 1150 Vec ll,rr; 1151 const PetscScalar *Y,*X; 1152 PetscScalar *x,*y; 1153 1154 PetscFunctionBegin; 1155 if (l) { 1156 ll = matis->y; 1157 PetscCall(VecGetArrayRead(l,&Y)); 1158 PetscCall(VecGetArray(ll,&y)); 1159 PetscCall(PetscSFBcastBegin(matis->sf,MPIU_SCALAR,Y,y,MPI_REPLACE)); 1160 } else { 1161 ll = NULL; 1162 } 1163 if (r) { 1164 rr = matis->x; 1165 PetscCall(VecGetArrayRead(r,&X)); 1166 PetscCall(VecGetArray(rr,&x)); 1167 PetscCall(PetscSFBcastBegin(matis->csf,MPIU_SCALAR,X,x,MPI_REPLACE)); 1168 } else { 1169 rr = NULL; 1170 } 1171 if (ll) { 1172 PetscCall(PetscSFBcastEnd(matis->sf,MPIU_SCALAR,Y,y,MPI_REPLACE)); 1173 PetscCall(VecRestoreArrayRead(l,&Y)); 1174 PetscCall(VecRestoreArray(ll,&y)); 1175 } 1176 if (rr) { 1177 PetscCall(PetscSFBcastEnd(matis->csf,MPIU_SCALAR,X,x,MPI_REPLACE)); 1178 PetscCall(VecRestoreArrayRead(r,&X)); 1179 PetscCall(VecRestoreArray(rr,&x)); 1180 } 1181 PetscCall(MatDiagonalScale(matis->A,ll,rr)); 1182 PetscFunctionReturn(0); 1183 } 1184 1185 static PetscErrorCode MatGetInfo_IS(Mat A,MatInfoType flag,MatInfo *ginfo) 1186 { 1187 Mat_IS *matis = (Mat_IS*)A->data; 1188 MatInfo info; 1189 PetscLogDouble isend[6],irecv[6]; 1190 PetscInt bs; 1191 1192 PetscFunctionBegin; 1193 PetscCall(MatGetBlockSize(A,&bs)); 1194 if (matis->A->ops->getinfo) { 1195 PetscCall(MatGetInfo(matis->A,MAT_LOCAL,&info)); 1196 isend[0] = info.nz_used; 1197 isend[1] = info.nz_allocated; 1198 isend[2] = info.nz_unneeded; 1199 isend[3] = info.memory; 1200 isend[4] = info.mallocs; 1201 } else { 1202 isend[0] = 0.; 1203 isend[1] = 0.; 1204 isend[2] = 0.; 1205 isend[3] = 0.; 1206 isend[4] = 0.; 1207 } 1208 isend[5] = matis->A->num_ass; 1209 if (flag == MAT_LOCAL) { 1210 ginfo->nz_used = isend[0]; 1211 ginfo->nz_allocated = isend[1]; 1212 ginfo->nz_unneeded = isend[2]; 1213 ginfo->memory = isend[3]; 1214 ginfo->mallocs = isend[4]; 1215 ginfo->assemblies = isend[5]; 1216 } else if (flag == MAT_GLOBAL_MAX) { 1217 PetscCall(MPIU_Allreduce(isend,irecv,6,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)A))); 1218 1219 ginfo->nz_used = irecv[0]; 1220 ginfo->nz_allocated = irecv[1]; 1221 ginfo->nz_unneeded = irecv[2]; 1222 ginfo->memory = irecv[3]; 1223 ginfo->mallocs = irecv[4]; 1224 ginfo->assemblies = irecv[5]; 1225 } else if (flag == MAT_GLOBAL_SUM) { 1226 PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)A))); 1227 1228 ginfo->nz_used = irecv[0]; 1229 ginfo->nz_allocated = irecv[1]; 1230 ginfo->nz_unneeded = irecv[2]; 1231 ginfo->memory = irecv[3]; 1232 ginfo->mallocs = irecv[4]; 1233 ginfo->assemblies = A->num_ass; 1234 } 1235 ginfo->block_size = bs; 1236 ginfo->fill_ratio_given = 0; 1237 ginfo->fill_ratio_needed = 0; 1238 ginfo->factor_mallocs = 0; 1239 PetscFunctionReturn(0); 1240 } 1241 1242 static PetscErrorCode MatTranspose_IS(Mat A,MatReuse reuse,Mat *B) 1243 { 1244 Mat C,lC,lA; 1245 1246 PetscFunctionBegin; 1247 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 1248 ISLocalToGlobalMapping rl2g,cl2g; 1249 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 1250 PetscCall(MatSetSizes(C,A->cmap->n,A->rmap->n,A->cmap->N,A->rmap->N)); 1251 PetscCall(MatSetBlockSizes(C,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs))); 1252 PetscCall(MatSetType(C,MATIS)); 1253 PetscCall(MatGetLocalToGlobalMapping(A,&rl2g,&cl2g)); 1254 PetscCall(MatSetLocalToGlobalMapping(C,cl2g,rl2g)); 1255 } else C = *B; 1256 1257 /* perform local transposition */ 1258 PetscCall(MatISGetLocalMat(A,&lA)); 1259 PetscCall(MatTranspose(lA,MAT_INITIAL_MATRIX,&lC)); 1260 PetscCall(MatSetLocalToGlobalMapping(lC,lA->cmap->mapping,lA->rmap->mapping)); 1261 PetscCall(MatISSetLocalMat(C,lC)); 1262 PetscCall(MatDestroy(&lC)); 1263 1264 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 1265 *B = C; 1266 } else { 1267 PetscCall(MatHeaderMerge(A,&C)); 1268 } 1269 PetscCall(MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY)); 1270 PetscCall(MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY)); 1271 PetscFunctionReturn(0); 1272 } 1273 1274 static PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 1275 { 1276 Mat_IS *is = (Mat_IS*)A->data; 1277 1278 PetscFunctionBegin; 1279 if (D) { /* MatShift_IS pass D = NULL */ 1280 PetscCall(VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD)); 1281 PetscCall(VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD)); 1282 } 1283 PetscCall(VecPointwiseDivide(is->y,is->y,is->counter)); 1284 PetscCall(MatDiagonalSet(is->A,is->y,insmode)); 1285 PetscFunctionReturn(0); 1286 } 1287 1288 static PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 1289 { 1290 Mat_IS *is = (Mat_IS*)A->data; 1291 1292 PetscFunctionBegin; 1293 PetscCall(VecSet(is->y,a)); 1294 PetscCall(MatDiagonalSet_IS(A,NULL,ADD_VALUES)); 1295 PetscFunctionReturn(0); 1296 } 1297 1298 static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1299 { 1300 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 1301 1302 PetscFunctionBegin; 1303 PetscCheckFalse(m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %d: they are %" PetscInt_FMT " %" PetscInt_FMT,MATIS_MAX_ENTRIES_INSERTION,m,n); 1304 PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l)); 1305 PetscCall(ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l)); 1306 PetscCall(MatSetValuesLocal_IS(A,m,rows_l,n,cols_l,values,addv)); 1307 PetscFunctionReturn(0); 1308 } 1309 1310 static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1311 { 1312 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 1313 1314 PetscFunctionBegin; 1315 PetscCheckFalse(m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %d: they are %" PetscInt_FMT " %" PetscInt_FMT,MATIS_MAX_ENTRIES_INSERTION,m,n); 1316 PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l)); 1317 PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l)); 1318 PetscCall(MatSetValuesBlockedLocal_IS(A,m,rows_l,n,cols_l,values,addv)); 1319 PetscFunctionReturn(0); 1320 } 1321 1322 static PetscErrorCode MatCreateSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 1323 { 1324 Mat locmat,newlocmat; 1325 Mat_IS *newmatis; 1326 const PetscInt *idxs; 1327 PetscInt i,m,n; 1328 1329 PetscFunctionBegin; 1330 if (scall == MAT_REUSE_MATRIX) { 1331 PetscBool ismatis; 1332 1333 PetscCall(PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis)); 1334 PetscCheck(ismatis,PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 1335 newmatis = (Mat_IS*)(*newmat)->data; 1336 PetscCheck(newmatis->getsub_ris,PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 1337 PetscCheck(newmatis->getsub_cis,PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 1338 } 1339 /* irow and icol may not have duplicate entries */ 1340 if (PetscDefined(USE_DEBUG)) { 1341 Vec rtest,ltest; 1342 const PetscScalar *array; 1343 1344 PetscCall(MatCreateVecs(mat,<est,&rtest)); 1345 PetscCall(ISGetLocalSize(irow,&n)); 1346 PetscCall(ISGetIndices(irow,&idxs)); 1347 for (i=0;i<n;i++) { 1348 PetscCall(VecSetValue(rtest,idxs[i],1.0,ADD_VALUES)); 1349 } 1350 PetscCall(VecAssemblyBegin(rtest)); 1351 PetscCall(VecAssemblyEnd(rtest)); 1352 PetscCall(VecGetLocalSize(rtest,&n)); 1353 PetscCall(VecGetOwnershipRange(rtest,&m,NULL)); 1354 PetscCall(VecGetArrayRead(rtest,&array)); 1355 for (i=0;i<n;i++) PetscCheckFalse(array[i] != 0. && array[i] != 1.,PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %" PetscInt_FMT " counted %" PetscInt_FMT " times! Irow may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i])); 1356 PetscCall(VecRestoreArrayRead(rtest,&array)); 1357 PetscCall(ISRestoreIndices(irow,&idxs)); 1358 PetscCall(ISGetLocalSize(icol,&n)); 1359 PetscCall(ISGetIndices(icol,&idxs)); 1360 for (i=0;i<n;i++) { 1361 PetscCall(VecSetValue(ltest,idxs[i],1.0,ADD_VALUES)); 1362 } 1363 PetscCall(VecAssemblyBegin(ltest)); 1364 PetscCall(VecAssemblyEnd(ltest)); 1365 PetscCall(VecGetLocalSize(ltest,&n)); 1366 PetscCall(VecGetOwnershipRange(ltest,&m,NULL)); 1367 PetscCall(VecGetArrayRead(ltest,&array)); 1368 for (i=0;i<n;i++) PetscCheckFalse(array[i] != 0. && array[i] != 1.,PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %" PetscInt_FMT " counted %" PetscInt_FMT " times! Icol may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i])); 1369 PetscCall(VecRestoreArrayRead(ltest,&array)); 1370 PetscCall(ISRestoreIndices(icol,&idxs)); 1371 PetscCall(VecDestroy(&rtest)); 1372 PetscCall(VecDestroy(<est)); 1373 } 1374 if (scall == MAT_INITIAL_MATRIX) { 1375 Mat_IS *matis = (Mat_IS*)mat->data; 1376 ISLocalToGlobalMapping rl2g; 1377 IS is; 1378 PetscInt *lidxs,*lgidxs,*newgidxs; 1379 PetscInt ll,newloc,irbs,icbs,arbs,acbs,rbs,cbs; 1380 PetscBool cong; 1381 MPI_Comm comm; 1382 1383 PetscCall(PetscObjectGetComm((PetscObject)mat,&comm)); 1384 PetscCall(MatGetBlockSizes(mat,&arbs,&acbs)); 1385 PetscCall(ISGetBlockSize(irow,&irbs)); 1386 PetscCall(ISGetBlockSize(icol,&icbs)); 1387 rbs = arbs == irbs ? irbs : 1; 1388 cbs = acbs == icbs ? icbs : 1; 1389 PetscCall(ISGetLocalSize(irow,&m)); 1390 PetscCall(ISGetLocalSize(icol,&n)); 1391 PetscCall(MatCreate(comm,newmat)); 1392 PetscCall(MatSetType(*newmat,MATIS)); 1393 PetscCall(MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE)); 1394 PetscCall(MatSetBlockSizes(*newmat,rbs,cbs)); 1395 /* communicate irow to their owners in the layout */ 1396 PetscCall(ISGetIndices(irow,&idxs)); 1397 PetscCall(PetscLayoutMapLocal(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs)); 1398 PetscCall(ISRestoreIndices(irow,&idxs)); 1399 PetscCall(PetscArrayzero(matis->sf_rootdata,matis->sf->nroots)); 1400 for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 1401 PetscCall(PetscFree(lidxs)); 1402 PetscCall(PetscFree(lgidxs)); 1403 PetscCall(PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata,MPI_REPLACE)); 1404 PetscCall(PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata,MPI_REPLACE)); 1405 for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 1406 PetscCall(PetscMalloc1(newloc,&newgidxs)); 1407 PetscCall(PetscMalloc1(newloc,&lidxs)); 1408 for (i=0,newloc=0;i<matis->sf->nleaves;i++) 1409 if (matis->sf_leafdata[i]) { 1410 lidxs[newloc] = i; 1411 newgidxs[newloc++] = matis->sf_leafdata[i]-1; 1412 } 1413 PetscCall(ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is)); 1414 PetscCall(ISLocalToGlobalMappingCreateIS(is,&rl2g)); 1415 PetscCall(ISLocalToGlobalMappingSetBlockSize(rl2g,rbs)); 1416 PetscCall(ISDestroy(&is)); 1417 /* local is to extract local submatrix */ 1418 newmatis = (Mat_IS*)(*newmat)->data; 1419 PetscCall(ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris)); 1420 PetscCall(MatHasCongruentLayouts(mat,&cong)); 1421 if (cong && irow == icol && matis->csf == matis->sf) { 1422 PetscCall(MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g)); 1423 PetscCall(PetscObjectReference((PetscObject)newmatis->getsub_ris)); 1424 newmatis->getsub_cis = newmatis->getsub_ris; 1425 } else { 1426 ISLocalToGlobalMapping cl2g; 1427 1428 /* communicate icol to their owners in the layout */ 1429 PetscCall(ISGetIndices(icol,&idxs)); 1430 PetscCall(PetscLayoutMapLocal(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs)); 1431 PetscCall(ISRestoreIndices(icol,&idxs)); 1432 PetscCall(PetscArrayzero(matis->csf_rootdata,matis->csf->nroots)); 1433 for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 1434 PetscCall(PetscFree(lidxs)); 1435 PetscCall(PetscFree(lgidxs)); 1436 PetscCall(PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata,MPI_REPLACE)); 1437 PetscCall(PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata,MPI_REPLACE)); 1438 for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 1439 PetscCall(PetscMalloc1(newloc,&newgidxs)); 1440 PetscCall(PetscMalloc1(newloc,&lidxs)); 1441 for (i=0,newloc=0;i<matis->csf->nleaves;i++) 1442 if (matis->csf_leafdata[i]) { 1443 lidxs[newloc] = i; 1444 newgidxs[newloc++] = matis->csf_leafdata[i]-1; 1445 } 1446 PetscCall(ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is)); 1447 PetscCall(ISLocalToGlobalMappingCreateIS(is,&cl2g)); 1448 PetscCall(ISLocalToGlobalMappingSetBlockSize(cl2g,cbs)); 1449 PetscCall(ISDestroy(&is)); 1450 /* local is to extract local submatrix */ 1451 PetscCall(ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis)); 1452 PetscCall(MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g)); 1453 PetscCall(ISLocalToGlobalMappingDestroy(&cl2g)); 1454 } 1455 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 1456 } else { 1457 PetscCall(MatISGetLocalMat(*newmat,&newlocmat)); 1458 } 1459 PetscCall(MatISGetLocalMat(mat,&locmat)); 1460 newmatis = (Mat_IS*)(*newmat)->data; 1461 PetscCall(MatCreateSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat)); 1462 if (scall == MAT_INITIAL_MATRIX) { 1463 PetscCall(MatISSetLocalMat(*newmat,newlocmat)); 1464 PetscCall(MatDestroy(&newlocmat)); 1465 } 1466 PetscCall(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY)); 1467 PetscCall(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY)); 1468 PetscFunctionReturn(0); 1469 } 1470 1471 static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 1472 { 1473 Mat_IS *a = (Mat_IS*)A->data,*b; 1474 PetscBool ismatis; 1475 1476 PetscFunctionBegin; 1477 PetscCall(PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis)); 1478 PetscCheck(ismatis,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 1479 b = (Mat_IS*)B->data; 1480 PetscCall(MatCopy(a->A,b->A,str)); 1481 PetscCall(PetscObjectStateIncrease((PetscObject)B)); 1482 PetscFunctionReturn(0); 1483 } 1484 1485 static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 1486 { 1487 Vec v; 1488 const PetscScalar *array; 1489 PetscInt i,n; 1490 1491 PetscFunctionBegin; 1492 *missing = PETSC_FALSE; 1493 PetscCall(MatCreateVecs(A,NULL,&v)); 1494 PetscCall(MatGetDiagonal(A,v)); 1495 PetscCall(VecGetLocalSize(v,&n)); 1496 PetscCall(VecGetArrayRead(v,&array)); 1497 for (i=0;i<n;i++) if (array[i] == 0.) break; 1498 PetscCall(VecRestoreArrayRead(v,&array)); 1499 PetscCall(VecDestroy(&v)); 1500 if (i != n) *missing = PETSC_TRUE; 1501 if (d) { 1502 *d = -1; 1503 if (*missing) { 1504 PetscInt rstart; 1505 PetscCall(MatGetOwnershipRange(A,&rstart,NULL)); 1506 *d = i+rstart; 1507 } 1508 } 1509 PetscFunctionReturn(0); 1510 } 1511 1512 static PetscErrorCode MatISSetUpSF_IS(Mat B) 1513 { 1514 Mat_IS *matis = (Mat_IS*)(B->data); 1515 const PetscInt *gidxs; 1516 PetscInt nleaves; 1517 1518 PetscFunctionBegin; 1519 if (matis->sf) PetscFunctionReturn(0); 1520 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf)); 1521 PetscCall(ISLocalToGlobalMappingGetIndices(matis->rmapping,&gidxs)); 1522 PetscCall(ISLocalToGlobalMappingGetSize(matis->rmapping,&nleaves)); 1523 PetscCall(PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs)); 1524 PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->rmapping,&gidxs)); 1525 PetscCall(PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata)); 1526 if (matis->rmapping != matis->cmapping) { /* setup SF for columns */ 1527 PetscCall(ISLocalToGlobalMappingGetSize(matis->cmapping,&nleaves)); 1528 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf)); 1529 PetscCall(ISLocalToGlobalMappingGetIndices(matis->cmapping,&gidxs)); 1530 PetscCall(PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs)); 1531 PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->cmapping,&gidxs)); 1532 PetscCall(PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata)); 1533 } else { 1534 matis->csf = matis->sf; 1535 matis->csf_leafdata = matis->sf_leafdata; 1536 matis->csf_rootdata = matis->sf_rootdata; 1537 } 1538 PetscFunctionReturn(0); 1539 } 1540 1541 /*@ 1542 MatISStoreL2L - Store local-to-local operators during the Galerkin process of MatPtAP. 1543 1544 Collective 1545 1546 Input Parameters: 1547 + A - the matrix 1548 - store - the boolean flag 1549 1550 Level: advanced 1551 1552 Notes: 1553 1554 .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatPtAP() 1555 @*/ 1556 PetscErrorCode MatISStoreL2L(Mat A, PetscBool store) 1557 { 1558 PetscFunctionBegin; 1559 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1560 PetscValidType(A,1); 1561 PetscValidLogicalCollectiveBool(A,store,2); 1562 PetscTryMethod(A,"MatISStoreL2L_C",(Mat,PetscBool),(A,store)); 1563 PetscFunctionReturn(0); 1564 } 1565 1566 static PetscErrorCode MatISStoreL2L_IS(Mat A, PetscBool store) 1567 { 1568 Mat_IS *matis = (Mat_IS*)(A->data); 1569 1570 PetscFunctionBegin; 1571 matis->storel2l = store; 1572 if (!store) { 1573 PetscCall(PetscObjectCompose((PetscObject)(A),"_MatIS_PtAP_l2l",NULL)); 1574 } 1575 PetscFunctionReturn(0); 1576 } 1577 1578 /*@ 1579 MatISFixLocalEmpty - Compress out zero local rows from the local matrices 1580 1581 Collective 1582 1583 Input Parameters: 1584 + A - the matrix 1585 - fix - the boolean flag 1586 1587 Level: advanced 1588 1589 Notes: When fix is true, new local matrices and l2g maps are generated during the final assembly process. 1590 1591 .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatAssemblyEnd(), MAT_FINAL_ASSEMBLY 1592 @*/ 1593 PetscErrorCode MatISFixLocalEmpty(Mat A, PetscBool fix) 1594 { 1595 PetscFunctionBegin; 1596 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1597 PetscValidType(A,1); 1598 PetscValidLogicalCollectiveBool(A,fix,2); 1599 PetscTryMethod(A,"MatISFixLocalEmpty_C",(Mat,PetscBool),(A,fix)); 1600 PetscFunctionReturn(0); 1601 } 1602 1603 static PetscErrorCode MatISFixLocalEmpty_IS(Mat A, PetscBool fix) 1604 { 1605 Mat_IS *matis = (Mat_IS*)(A->data); 1606 1607 PetscFunctionBegin; 1608 matis->locempty = fix; 1609 PetscFunctionReturn(0); 1610 } 1611 1612 /*@ 1613 MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 1614 1615 Collective 1616 1617 Input Parameters: 1618 + B - the matrix 1619 . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 1620 (same value is used for all local rows) 1621 . d_nnz - array containing the number of nonzeros in the various rows of the 1622 DIAGONAL portion of the local submatrix (possibly different for each row) 1623 or NULL, if d_nz is used to specify the nonzero structure. 1624 The size of this array is equal to the number of local rows, i.e 'm'. 1625 For matrices that will be factored, you must leave room for (and set) 1626 the diagonal entry even if it is zero. 1627 . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 1628 submatrix (same value is used for all local rows). 1629 - o_nnz - array containing the number of nonzeros in the various rows of the 1630 OFF-DIAGONAL portion of the local submatrix (possibly different for 1631 each row) or NULL, if o_nz is used to specify the nonzero 1632 structure. The size of this array is equal to the number 1633 of local rows, i.e 'm'. 1634 1635 If the *_nnz parameter is given then the *_nz parameter is ignored 1636 1637 Level: intermediate 1638 1639 Notes: 1640 This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 1641 from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 1642 matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 1643 1644 .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 1645 @*/ 1646 PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1647 { 1648 PetscFunctionBegin; 1649 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 1650 PetscValidType(B,1); 1651 PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz)); 1652 PetscFunctionReturn(0); 1653 } 1654 1655 /* this is used by DMDA */ 1656 PETSC_EXTERN PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1657 { 1658 Mat_IS *matis = (Mat_IS*)(B->data); 1659 PetscInt bs,i,nlocalcols; 1660 1661 PetscFunctionBegin; 1662 PetscCall(MatSetUp(B)); 1663 if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz; 1664 else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 1665 1666 if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz; 1667 else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 1668 1669 PetscCall(PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata,MPI_REPLACE)); 1670 PetscCall(MatGetSize(matis->A,NULL,&nlocalcols)); 1671 PetscCall(MatGetBlockSize(matis->A,&bs)); 1672 PetscCall(PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata,MPI_REPLACE)); 1673 1674 for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 1675 PetscCall(MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata)); 1676 #if defined(PETSC_HAVE_HYPRE) 1677 PetscCall(MatHYPRESetPreallocation(matis->A,0,matis->sf_leafdata,0,NULL)); 1678 #endif 1679 1680 for (i=0;i<matis->sf->nleaves/bs;i++) { 1681 PetscInt b; 1682 1683 matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 1684 for (b=1;b<bs;b++) { 1685 matis->sf_leafdata[i] = PetscMax(matis->sf_leafdata[i],matis->sf_leafdata[i*bs+b]/bs); 1686 } 1687 } 1688 PetscCall(MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata)); 1689 1690 nlocalcols /= bs; 1691 for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols - i); 1692 PetscCall(MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata)); 1693 1694 /* for other matrix types */ 1695 PetscCall(MatSetUp(matis->A)); 1696 PetscFunctionReturn(0); 1697 } 1698 1699 PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 1700 { 1701 Mat_IS *matis = (Mat_IS*)(A->data); 1702 PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 1703 const PetscInt *global_indices_r,*global_indices_c; 1704 PetscInt i,j,bs,rows,cols; 1705 PetscInt lrows,lcols; 1706 PetscInt local_rows,local_cols; 1707 PetscMPIInt size; 1708 PetscBool isdense,issbaij; 1709 1710 PetscFunctionBegin; 1711 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 1712 PetscCall(MatGetSize(A,&rows,&cols)); 1713 PetscCall(MatGetBlockSize(A,&bs)); 1714 PetscCall(MatGetSize(matis->A,&local_rows,&local_cols)); 1715 PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense)); 1716 PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij)); 1717 PetscCall(ISLocalToGlobalMappingGetIndices(matis->rmapping,&global_indices_r)); 1718 if (matis->rmapping != matis->cmapping) { 1719 PetscCall(ISLocalToGlobalMappingGetIndices(matis->cmapping,&global_indices_c)); 1720 } else global_indices_c = global_indices_r; 1721 1722 if (issbaij) PetscCall(MatGetRowUpperTriangular(matis->A)); 1723 /* 1724 An SF reduce is needed to sum up properly on shared rows. 1725 Note that generally preallocation is not exact, since it overestimates nonzeros 1726 */ 1727 PetscCall(MatGetLocalSize(A,&lrows,&lcols)); 1728 MatPreallocateBegin(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz); 1729 /* All processes need to compute entire row ownership */ 1730 PetscCall(PetscMalloc1(rows,&row_ownership)); 1731 PetscCall(MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges)); 1732 for (i=0;i<size;i++) { 1733 for (j=mat_ranges[i];j<mat_ranges[i+1];j++) row_ownership[j] = i; 1734 } 1735 PetscCall(MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges)); 1736 1737 /* 1738 my_dnz and my_onz contains exact contribution to preallocation from each local mat 1739 then, they will be summed up properly. This way, preallocation is always sufficient 1740 */ 1741 PetscCall(PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz)); 1742 /* preallocation as a MATAIJ */ 1743 if (isdense) { /* special case for dense local matrices */ 1744 for (i=0;i<local_rows;i++) { 1745 PetscInt owner = row_ownership[global_indices_r[i]]; 1746 for (j=0;j<local_cols;j++) { 1747 PetscInt index_col = global_indices_c[j]; 1748 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1]) { /* diag block */ 1749 my_dnz[i] += 1; 1750 } else { /* offdiag block */ 1751 my_onz[i] += 1; 1752 } 1753 } 1754 } 1755 } else if (matis->A->ops->getrowij) { 1756 const PetscInt *ii,*jj,*jptr; 1757 PetscBool done; 1758 PetscCall(MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done)); 1759 PetscCheck(done,PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ"); 1760 jptr = jj; 1761 for (i=0;i<local_rows;i++) { 1762 PetscInt index_row = global_indices_r[i]; 1763 for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 1764 PetscInt owner = row_ownership[index_row]; 1765 PetscInt index_col = global_indices_c[*jptr]; 1766 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1]) { /* diag block */ 1767 my_dnz[i] += 1; 1768 } else { /* offdiag block */ 1769 my_onz[i] += 1; 1770 } 1771 /* same as before, interchanging rows and cols */ 1772 if (issbaij && index_col != index_row) { 1773 owner = row_ownership[index_col]; 1774 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1]) { 1775 my_dnz[*jptr] += 1; 1776 } else { 1777 my_onz[*jptr] += 1; 1778 } 1779 } 1780 } 1781 } 1782 PetscCall(MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done)); 1783 PetscCheck(done,PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ"); 1784 } else { /* loop over rows and use MatGetRow */ 1785 for (i=0;i<local_rows;i++) { 1786 const PetscInt *cols; 1787 PetscInt ncols,index_row = global_indices_r[i]; 1788 PetscCall(MatGetRow(matis->A,i,&ncols,&cols,NULL)); 1789 for (j=0;j<ncols;j++) { 1790 PetscInt owner = row_ownership[index_row]; 1791 PetscInt index_col = global_indices_c[cols[j]]; 1792 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1]) { /* diag block */ 1793 my_dnz[i] += 1; 1794 } else { /* offdiag block */ 1795 my_onz[i] += 1; 1796 } 1797 /* same as before, interchanging rows and cols */ 1798 if (issbaij && index_col != index_row) { 1799 owner = row_ownership[index_col]; 1800 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1]) { 1801 my_dnz[cols[j]] += 1; 1802 } else { 1803 my_onz[cols[j]] += 1; 1804 } 1805 } 1806 } 1807 PetscCall(MatRestoreRow(matis->A,i,&ncols,&cols,NULL)); 1808 } 1809 } 1810 if (global_indices_c != global_indices_r) { 1811 PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->cmapping,&global_indices_c)); 1812 } 1813 PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->rmapping,&global_indices_r)); 1814 PetscCall(PetscFree(row_ownership)); 1815 1816 /* Reduce my_dnz and my_onz */ 1817 if (maxreduce) { 1818 PetscCall(PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX)); 1819 PetscCall(PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX)); 1820 PetscCall(PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX)); 1821 PetscCall(PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX)); 1822 } else { 1823 PetscCall(PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM)); 1824 PetscCall(PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM)); 1825 PetscCall(PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM)); 1826 PetscCall(PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM)); 1827 } 1828 PetscCall(PetscFree2(my_dnz,my_onz)); 1829 1830 /* Resize preallocation if overestimated */ 1831 for (i=0;i<lrows;i++) { 1832 dnz[i] = PetscMin(dnz[i],lcols); 1833 onz[i] = PetscMin(onz[i],cols-lcols); 1834 } 1835 1836 /* Set preallocation */ 1837 PetscCall(MatSetBlockSizesFromMats(B,A,A)); 1838 PetscCall(MatSeqAIJSetPreallocation(B,0,dnz)); 1839 PetscCall(MatMPIAIJSetPreallocation(B,0,dnz,0,onz)); 1840 for (i=0;i<lrows;i+=bs) { 1841 PetscInt b, d = dnz[i],o = onz[i]; 1842 1843 for (b=1;b<bs;b++) { 1844 d = PetscMax(d,dnz[i+b]); 1845 o = PetscMax(o,onz[i+b]); 1846 } 1847 dnz[i/bs] = PetscMin(d/bs + d%bs,lcols/bs); 1848 onz[i/bs] = PetscMin(o/bs + o%bs,(cols-lcols)/bs); 1849 } 1850 PetscCall(MatSeqBAIJSetPreallocation(B,bs,0,dnz)); 1851 PetscCall(MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz)); 1852 PetscCall(MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz)); 1853 MatPreallocateEnd(dnz,onz); 1854 if (issbaij) PetscCall(MatRestoreRowUpperTriangular(matis->A)); 1855 PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 1856 PetscFunctionReturn(0); 1857 } 1858 1859 PETSC_INTERN PetscErrorCode MatConvert_IS_XAIJ(Mat mat, MatType mtype, MatReuse reuse, Mat *M) 1860 { 1861 Mat_IS *matis = (Mat_IS*)(mat->data); 1862 Mat local_mat,MT; 1863 PetscInt rbs,cbs,rows,cols,lrows,lcols; 1864 PetscInt local_rows,local_cols; 1865 PetscBool isseqdense,isseqsbaij,isseqaij,isseqbaij; 1866 PetscMPIInt size; 1867 const PetscScalar *array; 1868 1869 PetscFunctionBegin; 1870 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size)); 1871 if (size == 1 && mat->rmap->N == matis->A->rmap->N && mat->cmap->N == matis->A->cmap->N) { 1872 Mat B; 1873 IS irows = NULL,icols = NULL; 1874 PetscInt rbs,cbs; 1875 1876 PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->rmapping,&rbs)); 1877 PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->cmapping,&cbs)); 1878 if (reuse != MAT_REUSE_MATRIX) { /* check if l2g maps are one-to-one */ 1879 IS rows,cols; 1880 const PetscInt *ridxs,*cidxs; 1881 PetscInt i,nw,*work; 1882 1883 PetscCall(ISLocalToGlobalMappingGetBlockIndices(matis->rmapping,&ridxs)); 1884 PetscCall(ISLocalToGlobalMappingGetSize(matis->rmapping,&nw)); 1885 nw = nw/rbs; 1886 PetscCall(PetscCalloc1(nw,&work)); 1887 for (i=0;i<nw;i++) work[ridxs[i]] += 1; 1888 for (i=0;i<nw;i++) if (!work[i] || work[i] > 1) break; 1889 if (i == nw) { 1890 PetscCall(ISCreateBlock(PETSC_COMM_SELF,rbs,nw,ridxs,PETSC_USE_POINTER,&rows)); 1891 PetscCall(ISSetPermutation(rows)); 1892 PetscCall(ISInvertPermutation(rows,PETSC_DECIDE,&irows)); 1893 PetscCall(ISDestroy(&rows)); 1894 } 1895 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(matis->rmapping,&ridxs)); 1896 PetscCall(PetscFree(work)); 1897 if (irows && matis->rmapping != matis->cmapping) { 1898 PetscCall(ISLocalToGlobalMappingGetBlockIndices(matis->cmapping,&cidxs)); 1899 PetscCall(ISLocalToGlobalMappingGetSize(matis->cmapping,&nw)); 1900 nw = nw/cbs; 1901 PetscCall(PetscCalloc1(nw,&work)); 1902 for (i=0;i<nw;i++) work[cidxs[i]] += 1; 1903 for (i=0;i<nw;i++) if (!work[i] || work[i] > 1) break; 1904 if (i == nw) { 1905 PetscCall(ISCreateBlock(PETSC_COMM_SELF,cbs,nw,cidxs,PETSC_USE_POINTER,&cols)); 1906 PetscCall(ISSetPermutation(cols)); 1907 PetscCall(ISInvertPermutation(cols,PETSC_DECIDE,&icols)); 1908 PetscCall(ISDestroy(&cols)); 1909 } 1910 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(matis->cmapping,&cidxs)); 1911 PetscCall(PetscFree(work)); 1912 } else if (irows) { 1913 PetscCall(PetscObjectReference((PetscObject)irows)); 1914 icols = irows; 1915 } 1916 } else { 1917 PetscCall(PetscObjectQuery((PetscObject)(*M),"_MatIS_IS_XAIJ_irows",(PetscObject*)&irows)); 1918 PetscCall(PetscObjectQuery((PetscObject)(*M),"_MatIS_IS_XAIJ_icols",(PetscObject*)&icols)); 1919 if (irows) PetscCall(PetscObjectReference((PetscObject)irows)); 1920 if (icols) PetscCall(PetscObjectReference((PetscObject)icols)); 1921 } 1922 if (!irows || !icols) { 1923 PetscCall(ISDestroy(&icols)); 1924 PetscCall(ISDestroy(&irows)); 1925 goto general_assembly; 1926 } 1927 PetscCall(MatConvert(matis->A,mtype,MAT_INITIAL_MATRIX,&B)); 1928 if (reuse != MAT_INPLACE_MATRIX) { 1929 PetscCall(MatCreateSubMatrix(B,irows,icols,reuse,M)); 1930 PetscCall(PetscObjectCompose((PetscObject)(*M),"_MatIS_IS_XAIJ_irows",(PetscObject)irows)); 1931 PetscCall(PetscObjectCompose((PetscObject)(*M),"_MatIS_IS_XAIJ_icols",(PetscObject)icols)); 1932 } else { 1933 Mat C; 1934 1935 PetscCall(MatCreateSubMatrix(B,irows,icols,MAT_INITIAL_MATRIX,&C)); 1936 PetscCall(MatHeaderReplace(mat,&C)); 1937 } 1938 PetscCall(MatDestroy(&B)); 1939 PetscCall(ISDestroy(&icols)); 1940 PetscCall(ISDestroy(&irows)); 1941 PetscFunctionReturn(0); 1942 } 1943 general_assembly: 1944 PetscCall(MatGetSize(mat,&rows,&cols)); 1945 PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->rmapping,&rbs)); 1946 PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->cmapping,&cbs)); 1947 PetscCall(MatGetLocalSize(mat,&lrows,&lcols)); 1948 PetscCall(MatGetSize(matis->A,&local_rows,&local_cols)); 1949 PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isseqdense)); 1950 PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij)); 1951 PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQBAIJ,&isseqbaij)); 1952 PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&isseqsbaij)); 1953 PetscCheckFalse(!isseqdense && !isseqaij && !isseqbaij && !isseqsbaij,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(matis->A))->type_name); 1954 if (PetscDefined (USE_DEBUG)) { 1955 PetscBool lb[4],bb[4]; 1956 1957 lb[0] = isseqdense; 1958 lb[1] = isseqaij; 1959 lb[2] = isseqbaij; 1960 lb[3] = isseqsbaij; 1961 PetscCall(MPIU_Allreduce(lb,bb,4,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat))); 1962 PetscCheckFalse(!bb[0] && !bb[1] && !bb[2] && !bb[3],PETSC_COMM_SELF,PETSC_ERR_SUP,"Local matrices must have the same type"); 1963 } 1964 1965 if (reuse != MAT_REUSE_MATRIX) { 1966 PetscCall(MatCreate(PetscObjectComm((PetscObject)mat),&MT)); 1967 PetscCall(MatSetSizes(MT,lrows,lcols,rows,cols)); 1968 PetscCall(MatSetType(MT,mtype)); 1969 PetscCall(MatSetBlockSizes(MT,rbs,cbs)); 1970 PetscCall(MatISSetMPIXAIJPreallocation_Private(mat,MT,PETSC_FALSE)); 1971 } else { 1972 PetscInt mrbs,mcbs,mrows,mcols,mlrows,mlcols; 1973 1974 /* some checks */ 1975 MT = *M; 1976 PetscCall(MatGetBlockSizes(MT,&mrbs,&mcbs)); 1977 PetscCall(MatGetSize(MT,&mrows,&mcols)); 1978 PetscCall(MatGetLocalSize(MT,&mlrows,&mlcols)); 1979 PetscCheck(mrows == rows,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%" PetscInt_FMT " != %" PetscInt_FMT ")",rows,mrows); 1980 PetscCheck(mcols == cols,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%" PetscInt_FMT " != %" PetscInt_FMT ")",cols,mcols); 1981 PetscCheck(mlrows == lrows,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%" PetscInt_FMT " != %" PetscInt_FMT ")",lrows,mlrows); 1982 PetscCheck(mlcols == lcols,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%" PetscInt_FMT " != %" PetscInt_FMT ")",lcols,mlcols); 1983 PetscCheck(mrbs == rbs,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong row block size (%" PetscInt_FMT " != %" PetscInt_FMT ")",rbs,mrbs); 1984 PetscCheck(mcbs == cbs,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong col block size (%" PetscInt_FMT " != %" PetscInt_FMT ")",cbs,mcbs); 1985 PetscCall(MatZeroEntries(MT)); 1986 } 1987 1988 if (isseqsbaij || isseqbaij) { 1989 PetscCall(MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&local_mat)); 1990 isseqaij = PETSC_TRUE; 1991 } else { 1992 PetscCall(PetscObjectReference((PetscObject)matis->A)); 1993 local_mat = matis->A; 1994 } 1995 1996 /* Set values */ 1997 PetscCall(MatSetLocalToGlobalMapping(MT,matis->rmapping,matis->cmapping)); 1998 if (isseqdense) { /* special case for dense local matrices */ 1999 PetscInt i,*dummy; 2000 2001 PetscCall(PetscMalloc1(PetscMax(local_rows,local_cols),&dummy)); 2002 for (i=0;i<PetscMax(local_rows,local_cols);i++) dummy[i] = i; 2003 PetscCall(MatSetOption(MT,MAT_ROW_ORIENTED,PETSC_FALSE)); 2004 PetscCall(MatDenseGetArrayRead(local_mat,&array)); 2005 PetscCall(MatSetValuesLocal(MT,local_rows,dummy,local_cols,dummy,array,ADD_VALUES)); 2006 PetscCall(MatDenseRestoreArrayRead(local_mat,&array)); 2007 PetscCall(PetscFree(dummy)); 2008 } else if (isseqaij) { 2009 const PetscInt *blocks; 2010 PetscInt i,nvtxs,*xadj,*adjncy, nb; 2011 PetscBool done; 2012 PetscScalar *sarray; 2013 2014 PetscCall(MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done)); 2015 PetscCheck(done,PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ"); 2016 PetscCall(MatSeqAIJGetArray(local_mat,&sarray)); 2017 PetscCall(MatGetVariableBlockSizes(local_mat,&nb,&blocks)); 2018 if (nb) { /* speed up assembly for special blocked matrices (used by BDDC) */ 2019 PetscInt sum; 2020 2021 for (i=0,sum=0;i<nb;i++) sum += blocks[i]; 2022 if (sum == nvtxs) { 2023 PetscInt r; 2024 2025 for (i=0,r=0;i<nb;i++) { 2026 PetscAssert(blocks[i] == xadj[r+1] - xadj[r],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid block sizes prescribed for block %" PetscInt_FMT ": expected %" PetscInt_FMT ", got %" PetscInt_FMT,i,blocks[i],xadj[r+1] - xadj[r]); 2027 PetscCall(MatSetValuesLocal(MT,blocks[i],adjncy+xadj[r],blocks[i],adjncy+xadj[r],sarray+xadj[r],ADD_VALUES)); 2028 r += blocks[i]; 2029 } 2030 } else { 2031 for (i=0;i<nvtxs;i++) { 2032 PetscCall(MatSetValuesLocal(MT,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],sarray+xadj[i],ADD_VALUES)); 2033 } 2034 } 2035 } else { 2036 for (i=0;i<nvtxs;i++) { 2037 PetscCall(MatSetValuesLocal(MT,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],sarray+xadj[i],ADD_VALUES)); 2038 } 2039 } 2040 PetscCall(MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done)); 2041 PetscCheck(done,PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ"); 2042 PetscCall(MatSeqAIJRestoreArray(local_mat,&sarray)); 2043 } else { /* very basic values insertion for all other matrix types */ 2044 PetscInt i; 2045 2046 for (i=0;i<local_rows;i++) { 2047 PetscInt j; 2048 const PetscInt *local_indices_cols; 2049 2050 PetscCall(MatGetRow(local_mat,i,&j,&local_indices_cols,&array)); 2051 PetscCall(MatSetValuesLocal(MT,1,&i,j,local_indices_cols,array,ADD_VALUES)); 2052 PetscCall(MatRestoreRow(local_mat,i,&j,&local_indices_cols,&array)); 2053 } 2054 } 2055 PetscCall(MatAssemblyBegin(MT,MAT_FINAL_ASSEMBLY)); 2056 PetscCall(MatDestroy(&local_mat)); 2057 PetscCall(MatAssemblyEnd(MT,MAT_FINAL_ASSEMBLY)); 2058 if (isseqdense) { 2059 PetscCall(MatSetOption(MT,MAT_ROW_ORIENTED,PETSC_TRUE)); 2060 } 2061 if (reuse == MAT_INPLACE_MATRIX) { 2062 PetscCall(MatHeaderReplace(mat,&MT)); 2063 } else if (reuse == MAT_INITIAL_MATRIX) { 2064 *M = MT; 2065 } 2066 PetscFunctionReturn(0); 2067 } 2068 2069 /*@ 2070 MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 2071 2072 Input Parameters: 2073 + mat - the matrix (should be of type MATIS) 2074 - reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 2075 2076 Output Parameter: 2077 . newmat - the matrix in AIJ format 2078 2079 Level: developer 2080 2081 Notes: 2082 This function has been deprecated and it will be removed in future releases. Update your code to use the MatConvert() interface. 2083 2084 .seealso: MATIS, MatConvert() 2085 @*/ 2086 PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 2087 { 2088 PetscFunctionBegin; 2089 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2090 PetscValidLogicalCollectiveEnum(mat,reuse,2); 2091 PetscValidPointer(newmat,3); 2092 if (reuse == MAT_REUSE_MATRIX) { 2093 PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 2094 PetscCheckSameComm(mat,1,*newmat,3); 2095 PetscCheck(mat != *newmat,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 2096 } 2097 PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatType,MatReuse,Mat*),(mat,MATAIJ,reuse,newmat)); 2098 PetscFunctionReturn(0); 2099 } 2100 2101 static PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 2102 { 2103 Mat_IS *matis = (Mat_IS*)(mat->data); 2104 PetscInt rbs,cbs,m,n,M,N; 2105 Mat B,localmat; 2106 2107 PetscFunctionBegin; 2108 PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->rmap->mapping,&rbs)); 2109 PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->cmap->mapping,&cbs)); 2110 PetscCall(MatGetSize(mat,&M,&N)); 2111 PetscCall(MatGetLocalSize(mat,&m,&n)); 2112 PetscCall(MatCreate(PetscObjectComm((PetscObject)mat),&B)); 2113 PetscCall(MatSetSizes(B,m,n,M,N)); 2114 PetscCall(MatSetBlockSize(B,rbs == cbs ? rbs : 1)); 2115 PetscCall(MatSetType(B,MATIS)); 2116 PetscCall(MatISSetLocalMatType(B,matis->lmattype)); 2117 PetscCall(MatSetLocalToGlobalMapping(B,mat->rmap->mapping,mat->cmap->mapping)); 2118 PetscCall(MatDuplicate(matis->A,op,&localmat)); 2119 PetscCall(MatSetLocalToGlobalMapping(localmat,matis->A->rmap->mapping,matis->A->cmap->mapping)); 2120 PetscCall(MatISSetLocalMat(B,localmat)); 2121 PetscCall(MatDestroy(&localmat)); 2122 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2123 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 2124 *newmat = B; 2125 PetscFunctionReturn(0); 2126 } 2127 2128 static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 2129 { 2130 Mat_IS *matis = (Mat_IS*)A->data; 2131 PetscBool local_sym; 2132 2133 PetscFunctionBegin; 2134 PetscCall(MatIsHermitian(matis->A,tol,&local_sym)); 2135 PetscCall(MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 2136 PetscFunctionReturn(0); 2137 } 2138 2139 static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 2140 { 2141 Mat_IS *matis = (Mat_IS*)A->data; 2142 PetscBool local_sym; 2143 2144 PetscFunctionBegin; 2145 if (matis->rmapping != matis->cmapping) { 2146 *flg = PETSC_FALSE; 2147 PetscFunctionReturn(0); 2148 } 2149 PetscCall(MatIsSymmetric(matis->A,tol,&local_sym)); 2150 PetscCall(MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 2151 PetscFunctionReturn(0); 2152 } 2153 2154 static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A,PetscBool *flg) 2155 { 2156 Mat_IS *matis = (Mat_IS*)A->data; 2157 PetscBool local_sym; 2158 2159 PetscFunctionBegin; 2160 if (matis->rmapping != matis->cmapping) { 2161 *flg = PETSC_FALSE; 2162 PetscFunctionReturn(0); 2163 } 2164 PetscCall(MatIsStructurallySymmetric(matis->A,&local_sym)); 2165 PetscCall(MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 2166 PetscFunctionReturn(0); 2167 } 2168 2169 static PetscErrorCode MatDestroy_IS(Mat A) 2170 { 2171 Mat_IS *b = (Mat_IS*)A->data; 2172 2173 PetscFunctionBegin; 2174 PetscCall(PetscFree(b->bdiag)); 2175 PetscCall(PetscFree(b->lmattype)); 2176 PetscCall(MatDestroy(&b->A)); 2177 PetscCall(VecScatterDestroy(&b->cctx)); 2178 PetscCall(VecScatterDestroy(&b->rctx)); 2179 PetscCall(VecDestroy(&b->x)); 2180 PetscCall(VecDestroy(&b->y)); 2181 PetscCall(VecDestroy(&b->counter)); 2182 PetscCall(ISDestroy(&b->getsub_ris)); 2183 PetscCall(ISDestroy(&b->getsub_cis)); 2184 if (b->sf != b->csf) { 2185 PetscCall(PetscSFDestroy(&b->csf)); 2186 PetscCall(PetscFree2(b->csf_rootdata,b->csf_leafdata)); 2187 } else b->csf = NULL; 2188 PetscCall(PetscSFDestroy(&b->sf)); 2189 PetscCall(PetscFree2(b->sf_rootdata,b->sf_leafdata)); 2190 PetscCall(ISLocalToGlobalMappingDestroy(&b->rmapping)); 2191 PetscCall(ISLocalToGlobalMappingDestroy(&b->cmapping)); 2192 PetscCall(PetscFree(A->data)); 2193 PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL)); 2194 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMatType_C",NULL)); 2195 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL)); 2196 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL)); 2197 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL)); 2198 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL)); 2199 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISStoreL2L_C",NULL)); 2200 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISFixLocalEmpty_C",NULL)); 2201 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalToGlobalMapping_C",NULL)); 2202 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_mpiaij_C",NULL)); 2203 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_mpibaij_C",NULL)); 2204 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_mpisbaij_C",NULL)); 2205 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_seqaij_C",NULL)); 2206 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_seqbaij_C",NULL)); 2207 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_seqsbaij_C",NULL)); 2208 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_aij_C",NULL)); 2209 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",NULL)); 2210 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetPreallocationCOO_C",NULL)); 2211 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetValuesCOO_C",NULL)); 2212 PetscFunctionReturn(0); 2213 } 2214 2215 static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 2216 { 2217 Mat_IS *is = (Mat_IS*)A->data; 2218 PetscScalar zero = 0.0; 2219 2220 PetscFunctionBegin; 2221 /* scatter the global vector x into the local work vector */ 2222 PetscCall(VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD)); 2223 PetscCall(VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD)); 2224 2225 /* multiply the local matrix */ 2226 PetscCall(MatMult(is->A,is->x,is->y)); 2227 2228 /* scatter product back into global memory */ 2229 PetscCall(VecSet(y,zero)); 2230 PetscCall(VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE)); 2231 PetscCall(VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE)); 2232 PetscFunctionReturn(0); 2233 } 2234 2235 static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 2236 { 2237 Vec temp_vec; 2238 2239 PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 2240 if (v3 != v2) { 2241 PetscCall(MatMult(A,v1,v3)); 2242 PetscCall(VecAXPY(v3,1.0,v2)); 2243 } else { 2244 PetscCall(VecDuplicate(v2,&temp_vec)); 2245 PetscCall(MatMult(A,v1,temp_vec)); 2246 PetscCall(VecAXPY(temp_vec,1.0,v2)); 2247 PetscCall(VecCopy(temp_vec,v3)); 2248 PetscCall(VecDestroy(&temp_vec)); 2249 } 2250 PetscFunctionReturn(0); 2251 } 2252 2253 static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 2254 { 2255 Mat_IS *is = (Mat_IS*)A->data; 2256 2257 PetscFunctionBegin; 2258 /* scatter the global vector x into the local work vector */ 2259 PetscCall(VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD)); 2260 PetscCall(VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD)); 2261 2262 /* multiply the local matrix */ 2263 PetscCall(MatMultTranspose(is->A,is->y,is->x)); 2264 2265 /* scatter product back into global vector */ 2266 PetscCall(VecSet(x,0)); 2267 PetscCall(VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE)); 2268 PetscCall(VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE)); 2269 PetscFunctionReturn(0); 2270 } 2271 2272 static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 2273 { 2274 Vec temp_vec; 2275 2276 PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 2277 if (v3 != v2) { 2278 PetscCall(MatMultTranspose(A,v1,v3)); 2279 PetscCall(VecAXPY(v3,1.0,v2)); 2280 } else { 2281 PetscCall(VecDuplicate(v2,&temp_vec)); 2282 PetscCall(MatMultTranspose(A,v1,temp_vec)); 2283 PetscCall(VecAXPY(temp_vec,1.0,v2)); 2284 PetscCall(VecCopy(temp_vec,v3)); 2285 PetscCall(VecDestroy(&temp_vec)); 2286 } 2287 PetscFunctionReturn(0); 2288 } 2289 2290 static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 2291 { 2292 Mat_IS *a = (Mat_IS*)A->data; 2293 PetscViewer sviewer; 2294 PetscBool isascii,view = PETSC_TRUE; 2295 2296 PetscFunctionBegin; 2297 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 2298 if (isascii) { 2299 PetscViewerFormat format; 2300 2301 PetscCall(PetscViewerGetFormat(viewer,&format)); 2302 if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE; 2303 } 2304 if (!view) PetscFunctionReturn(0); 2305 PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 2306 PetscCall(MatView(a->A,sviewer)); 2307 PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 2308 PetscCall(PetscViewerFlush(viewer)); 2309 PetscFunctionReturn(0); 2310 } 2311 2312 static PetscErrorCode MatInvertBlockDiagonal_IS(Mat mat,const PetscScalar **values) 2313 { 2314 Mat_IS *is = (Mat_IS*)mat->data; 2315 MPI_Datatype nodeType; 2316 const PetscScalar *lv; 2317 PetscInt bs; 2318 2319 PetscFunctionBegin; 2320 PetscCall(MatGetBlockSize(mat,&bs)); 2321 PetscCall(MatSetBlockSize(is->A,bs)); 2322 PetscCall(MatInvertBlockDiagonal(is->A,&lv)); 2323 if (!is->bdiag) { 2324 PetscCall(PetscMalloc1(bs*mat->rmap->n,&is->bdiag)); 2325 } 2326 PetscCallMPI(MPI_Type_contiguous(bs,MPIU_SCALAR,&nodeType)); 2327 PetscCallMPI(MPI_Type_commit(&nodeType)); 2328 PetscCall(PetscSFReduceBegin(is->sf,nodeType,lv,is->bdiag,MPI_REPLACE)); 2329 PetscCall(PetscSFReduceEnd(is->sf,nodeType,lv,is->bdiag,MPI_REPLACE)); 2330 PetscCallMPI(MPI_Type_free(&nodeType)); 2331 if (values) *values = is->bdiag; 2332 PetscFunctionReturn(0); 2333 } 2334 2335 static PetscErrorCode MatISSetUpScatters_Private(Mat A) 2336 { 2337 Vec cglobal,rglobal; 2338 IS from; 2339 Mat_IS *is = (Mat_IS*)A->data; 2340 PetscScalar sum; 2341 const PetscInt *garray; 2342 PetscInt nr,rbs,nc,cbs; 2343 VecType rtype; 2344 2345 PetscFunctionBegin; 2346 PetscCall(ISLocalToGlobalMappingGetSize(is->rmapping,&nr)); 2347 PetscCall(ISLocalToGlobalMappingGetBlockSize(is->rmapping,&rbs)); 2348 PetscCall(ISLocalToGlobalMappingGetSize(is->cmapping,&nc)); 2349 PetscCall(ISLocalToGlobalMappingGetBlockSize(is->cmapping,&cbs)); 2350 PetscCall(VecDestroy(&is->x)); 2351 PetscCall(VecDestroy(&is->y)); 2352 PetscCall(VecDestroy(&is->counter)); 2353 PetscCall(VecScatterDestroy(&is->rctx)); 2354 PetscCall(VecScatterDestroy(&is->cctx)); 2355 PetscCall(MatCreateVecs(is->A,&is->x,&is->y)); 2356 PetscCall(VecBindToCPU(is->y,PETSC_TRUE)); 2357 PetscCall(VecGetRootType_Private(is->y,&rtype)); 2358 PetscCall(PetscFree(A->defaultvectype)); 2359 PetscCall(PetscStrallocpy(rtype,&A->defaultvectype)); 2360 PetscCall(MatCreateVecs(A,&cglobal,&rglobal)); 2361 PetscCall(VecBindToCPU(rglobal,PETSC_TRUE)); 2362 PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->rmapping,&garray)); 2363 PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A),rbs,nr/rbs,garray,PETSC_USE_POINTER,&from)); 2364 PetscCall(VecScatterCreate(rglobal,from,is->y,NULL,&is->rctx)); 2365 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->rmapping,&garray)); 2366 PetscCall(ISDestroy(&from)); 2367 if (is->rmapping != is->cmapping) { 2368 PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->cmapping,&garray)); 2369 PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A),cbs,nc/cbs,garray,PETSC_USE_POINTER,&from)); 2370 PetscCall(VecScatterCreate(cglobal,from,is->x,NULL,&is->cctx)); 2371 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->cmapping,&garray)); 2372 PetscCall(ISDestroy(&from)); 2373 } else { 2374 PetscCall(PetscObjectReference((PetscObject)is->rctx)); 2375 is->cctx = is->rctx; 2376 } 2377 PetscCall(VecDestroy(&cglobal)); 2378 2379 /* interface counter vector (local) */ 2380 PetscCall(VecDuplicate(is->y,&is->counter)); 2381 PetscCall(VecBindToCPU(is->counter,PETSC_TRUE)); 2382 PetscCall(VecSet(is->y,1.)); 2383 PetscCall(VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE)); 2384 PetscCall(VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE)); 2385 PetscCall(VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD)); 2386 PetscCall(VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD)); 2387 PetscCall(VecBindToCPU(is->y,PETSC_FALSE)); 2388 PetscCall(VecBindToCPU(is->counter,PETSC_FALSE)); 2389 2390 /* special functions for block-diagonal matrices */ 2391 PetscCall(VecSum(rglobal,&sum)); 2392 A->ops->invertblockdiagonal = NULL; 2393 if ((PetscInt)(PetscRealPart(sum)) == A->rmap->N && A->rmap->N == A->cmap->N && is->rmapping == is->cmapping) A->ops->invertblockdiagonal = MatInvertBlockDiagonal_IS; 2394 PetscCall(VecDestroy(&rglobal)); 2395 2396 /* setup SF for general purpose shared indices based communications */ 2397 PetscCall(MatISSetUpSF_IS(A)); 2398 PetscFunctionReturn(0); 2399 } 2400 2401 static PetscErrorCode MatISFilterL2GMap(Mat A, ISLocalToGlobalMapping map, ISLocalToGlobalMapping *nmap, ISLocalToGlobalMapping *lmap) 2402 { 2403 IS is; 2404 ISLocalToGlobalMappingType l2gtype; 2405 const PetscInt *idxs; 2406 PetscHSetI ht; 2407 PetscInt *nidxs; 2408 PetscInt i,n,bs,c; 2409 PetscBool flg[] = {PETSC_FALSE,PETSC_FALSE}; 2410 2411 PetscFunctionBegin; 2412 PetscCall(ISLocalToGlobalMappingGetSize(map,&n)); 2413 PetscCall(ISLocalToGlobalMappingGetBlockSize(map,&bs)); 2414 PetscCall(ISLocalToGlobalMappingGetBlockIndices(map,&idxs)); 2415 PetscCall(PetscHSetICreate(&ht)); 2416 PetscCall(PetscMalloc1(n/bs,&nidxs)); 2417 for (i=0,c=0;i<n/bs;i++) { 2418 PetscBool missing; 2419 if (idxs[i] < 0) { flg[0] = PETSC_TRUE; continue; } 2420 PetscCall(PetscHSetIQueryAdd(ht,idxs[i],&missing)); 2421 if (!missing) flg[1] = PETSC_TRUE; 2422 else nidxs[c++] = idxs[i]; 2423 } 2424 PetscCall(PetscHSetIDestroy(&ht)); 2425 PetscCall(MPIU_Allreduce(MPI_IN_PLACE,flg,2,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)A))); 2426 if (!flg[0] && !flg[1]) { /* Entries are all non negative and unique */ 2427 *nmap = NULL; 2428 *lmap = NULL; 2429 PetscCall(PetscFree(nidxs)); 2430 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(map,&idxs)); 2431 PetscFunctionReturn(0); 2432 } 2433 2434 /* New l2g map without negative or repeated indices */ 2435 PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A),bs,c,nidxs,PETSC_USE_POINTER,&is)); 2436 PetscCall(ISLocalToGlobalMappingCreateIS(is,nmap)); 2437 PetscCall(ISDestroy(&is)); 2438 PetscCall(ISLocalToGlobalMappingGetType(map,&l2gtype)); 2439 PetscCall(ISLocalToGlobalMappingSetType(*nmap,l2gtype)); 2440 2441 /* New local l2g map for repeated indices */ 2442 PetscCall(ISGlobalToLocalMappingApplyBlock(*nmap,IS_GTOLM_MASK,n/bs,idxs,NULL,nidxs)); 2443 PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,n/bs,nidxs,PETSC_USE_POINTER,&is)); 2444 PetscCall(ISLocalToGlobalMappingCreateIS(is,lmap)); 2445 PetscCall(ISDestroy(&is)); 2446 2447 PetscCall(PetscFree(nidxs)); 2448 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(map,&idxs)); 2449 PetscFunctionReturn(0); 2450 } 2451 2452 static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 2453 { 2454 Mat_IS *is = (Mat_IS*)A->data; 2455 ISLocalToGlobalMapping localrmapping = NULL, localcmapping = NULL; 2456 PetscBool cong, freem[] = { PETSC_FALSE, PETSC_FALSE }; 2457 PetscInt nr,rbs,nc,cbs; 2458 2459 PetscFunctionBegin; 2460 if (rmapping) PetscCheckSameComm(A,1,rmapping,2); 2461 if (cmapping) PetscCheckSameComm(A,1,cmapping,3); 2462 2463 PetscCall(ISLocalToGlobalMappingDestroy(&is->rmapping)); 2464 PetscCall(ISLocalToGlobalMappingDestroy(&is->cmapping)); 2465 PetscCall(PetscLayoutSetUp(A->rmap)); 2466 PetscCall(PetscLayoutSetUp(A->cmap)); 2467 PetscCall(MatHasCongruentLayouts(A,&cong)); 2468 2469 /* If NULL, local space matches global space */ 2470 if (!rmapping) { 2471 IS is; 2472 2473 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A),A->rmap->N,0,1,&is)); 2474 PetscCall(ISLocalToGlobalMappingCreateIS(is,&rmapping)); 2475 if (A->rmap->bs > 0) PetscCall(ISLocalToGlobalMappingSetBlockSize(rmapping,A->rmap->bs)); 2476 PetscCall(ISDestroy(&is)); 2477 freem[0] = PETSC_TRUE; 2478 if (!cmapping && cong && A->rmap->bs == A->cmap->bs) cmapping = rmapping; 2479 } else if (!is->islocalref) { /* check if the l2g map has negative or repeated entries */ 2480 PetscCall(MatISFilterL2GMap(A,rmapping,&is->rmapping,&localrmapping)); 2481 if (rmapping == cmapping) { 2482 PetscCall(PetscObjectReference((PetscObject)is->rmapping)); 2483 is->cmapping = is->rmapping; 2484 PetscCall(PetscObjectReference((PetscObject)localrmapping)); 2485 localcmapping = localrmapping; 2486 } 2487 } 2488 if (!cmapping) { 2489 IS is; 2490 2491 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A),A->cmap->N,0,1,&is)); 2492 PetscCall(ISLocalToGlobalMappingCreateIS(is,&cmapping)); 2493 if (A->cmap->bs > 0) PetscCall(ISLocalToGlobalMappingSetBlockSize(cmapping,A->cmap->bs)); 2494 PetscCall(ISDestroy(&is)); 2495 freem[1] = PETSC_TRUE; 2496 } else if (cmapping != rmapping && !is->islocalref) { /* check if the l2g map has negative or repeated entries */ 2497 PetscCall(MatISFilterL2GMap(A,cmapping,&is->cmapping,&localcmapping)); 2498 } 2499 if (!is->rmapping) { 2500 PetscCall(PetscObjectReference((PetscObject)rmapping)); 2501 is->rmapping = rmapping; 2502 } 2503 if (!is->cmapping) { 2504 PetscCall(PetscObjectReference((PetscObject)cmapping)); 2505 is->cmapping = cmapping; 2506 } 2507 2508 /* Clean up */ 2509 PetscCall(MatDestroy(&is->A)); 2510 if (is->csf != is->sf) { 2511 PetscCall(PetscSFDestroy(&is->csf)); 2512 PetscCall(PetscFree2(is->csf_rootdata,is->csf_leafdata)); 2513 } else is->csf = NULL; 2514 PetscCall(PetscSFDestroy(&is->sf)); 2515 PetscCall(PetscFree2(is->sf_rootdata,is->sf_leafdata)); 2516 PetscCall(PetscFree(is->bdiag)); 2517 2518 /* check if the two mappings are actually the same for square matrices since MATIS has some optimization for this case 2519 (DOLFIN passes 2 different objects) */ 2520 PetscCall(ISLocalToGlobalMappingGetSize(is->rmapping,&nr)); 2521 PetscCall(ISLocalToGlobalMappingGetBlockSize(is->rmapping,&rbs)); 2522 PetscCall(ISLocalToGlobalMappingGetSize(is->cmapping,&nc)); 2523 PetscCall(ISLocalToGlobalMappingGetBlockSize(is->cmapping,&cbs)); 2524 if (is->rmapping != is->cmapping && cong) { 2525 PetscBool same = PETSC_FALSE; 2526 if (nr == nc && cbs == rbs) { 2527 const PetscInt *idxs1,*idxs2; 2528 2529 PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->rmapping,&idxs1)); 2530 PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->cmapping,&idxs2)); 2531 PetscCall(PetscArraycmp(idxs1,idxs2,nr/rbs,&same)); 2532 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->rmapping,&idxs1)); 2533 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->cmapping,&idxs2)); 2534 } 2535 PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&same,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 2536 if (same) { 2537 PetscCall(ISLocalToGlobalMappingDestroy(&is->cmapping)); 2538 PetscCall(PetscObjectReference((PetscObject)is->rmapping)); 2539 is->cmapping = is->rmapping; 2540 } 2541 } 2542 PetscCall(PetscLayoutSetBlockSize(A->rmap,rbs)); 2543 PetscCall(PetscLayoutSetBlockSize(A->cmap,cbs)); 2544 /* Pass the user defined maps to the layout */ 2545 PetscCall(PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping)); 2546 PetscCall(PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping)); 2547 if (freem[0]) PetscCall(ISLocalToGlobalMappingDestroy(&rmapping)); 2548 if (freem[1]) PetscCall(ISLocalToGlobalMappingDestroy(&cmapping)); 2549 2550 /* Create the local matrix A */ 2551 PetscCall(MatCreate(PETSC_COMM_SELF,&is->A)); 2552 PetscCall(MatSetType(is->A,is->lmattype)); 2553 PetscCall(MatSetSizes(is->A,nr,nc,nr,nc)); 2554 PetscCall(MatSetBlockSizes(is->A,rbs,cbs)); 2555 PetscCall(MatSetOptionsPrefix(is->A,"is_")); 2556 PetscCall(MatAppendOptionsPrefix(is->A,((PetscObject)A)->prefix)); 2557 PetscCall(PetscLayoutSetUp(is->A->rmap)); 2558 PetscCall(PetscLayoutSetUp(is->A->cmap)); 2559 PetscCall(MatSetLocalToGlobalMapping(is->A,localrmapping,localcmapping)); 2560 PetscCall(ISLocalToGlobalMappingDestroy(&localrmapping)); 2561 PetscCall(ISLocalToGlobalMappingDestroy(&localcmapping)); 2562 2563 /* setup scatters and local vectors for MatMult */ 2564 if (!is->islocalref) PetscCall(MatISSetUpScatters_Private(A)); 2565 A->preallocated = PETSC_TRUE; 2566 PetscFunctionReturn(0); 2567 } 2568 2569 static PetscErrorCode MatSetUp_IS(Mat A) 2570 { 2571 ISLocalToGlobalMapping rmap, cmap; 2572 2573 PetscFunctionBegin; 2574 PetscCall(MatGetLocalToGlobalMapping(A,&rmap,&cmap)); 2575 if (!rmap && !cmap) { 2576 PetscCall(MatSetLocalToGlobalMapping(A,NULL,NULL)); 2577 } 2578 PetscFunctionReturn(0); 2579 } 2580 2581 static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 2582 { 2583 Mat_IS *is = (Mat_IS*)mat->data; 2584 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 2585 2586 PetscFunctionBegin; 2587 PetscCall(ISGlobalToLocalMappingApply(is->rmapping,IS_GTOLM_MASK,m,rows,&m,rows_l)); 2588 if (m != n || rows != cols || is->cmapping != is->rmapping) { 2589 PetscCall(ISGlobalToLocalMappingApply(is->cmapping,IS_GTOLM_MASK,n,cols,&n,cols_l)); 2590 PetscCall(MatSetValues(is->A,m,rows_l,n,cols_l,values,addv)); 2591 } else { 2592 PetscCall(MatSetValues(is->A,m,rows_l,m,rows_l,values,addv)); 2593 } 2594 PetscFunctionReturn(0); 2595 } 2596 2597 static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 2598 { 2599 Mat_IS *is = (Mat_IS*)mat->data; 2600 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 2601 2602 PetscFunctionBegin; 2603 PetscCall(ISGlobalToLocalMappingApplyBlock(is->rmapping,IS_GTOLM_MASK,m,rows,&m,rows_l)); 2604 if (m != n || rows != cols || is->cmapping != is->rmapping) { 2605 PetscCall(ISGlobalToLocalMappingApply(is->cmapping,IS_GTOLM_MASK,n,cols,&n,cols_l)); 2606 PetscCall(MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv)); 2607 } else { 2608 PetscCall(MatSetValuesBlocked(is->A,m,rows_l,n,rows_l,values,addv)); 2609 } 2610 PetscFunctionReturn(0); 2611 } 2612 2613 static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 2614 { 2615 Mat_IS *is = (Mat_IS*)A->data; 2616 2617 PetscFunctionBegin; 2618 if (is->A->rmap->mapping || is->A->cmap->mapping) { 2619 PetscCall(MatSetValuesLocal(is->A,m,rows,n,cols,values,addv)); 2620 } else { 2621 PetscCall(MatSetValues(is->A,m,rows,n,cols,values,addv)); 2622 } 2623 PetscFunctionReturn(0); 2624 } 2625 2626 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 2627 { 2628 Mat_IS *is = (Mat_IS*)A->data; 2629 2630 PetscFunctionBegin; 2631 if (is->A->rmap->mapping || is->A->cmap->mapping) { 2632 PetscCall(MatSetValuesBlockedLocal(is->A,m,rows,n,cols,values,addv)); 2633 } else { 2634 PetscCall(MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv)); 2635 } 2636 PetscFunctionReturn(0); 2637 } 2638 2639 static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns) 2640 { 2641 Mat_IS *is = (Mat_IS*)A->data; 2642 2643 PetscFunctionBegin; 2644 if (!n) { 2645 is->pure_neumann = PETSC_TRUE; 2646 } else { 2647 PetscInt i; 2648 is->pure_neumann = PETSC_FALSE; 2649 2650 if (columns) { 2651 PetscCall(MatZeroRowsColumns(is->A,n,rows,diag,NULL,NULL)); 2652 } else { 2653 PetscCall(MatZeroRows(is->A,n,rows,diag,NULL,NULL)); 2654 } 2655 if (diag != 0.) { 2656 const PetscScalar *array; 2657 PetscCall(VecGetArrayRead(is->counter,&array)); 2658 for (i=0; i<n; i++) { 2659 PetscCall(MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES)); 2660 } 2661 PetscCall(VecRestoreArrayRead(is->counter,&array)); 2662 } 2663 PetscCall(MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY)); 2664 PetscCall(MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY)); 2665 } 2666 PetscFunctionReturn(0); 2667 } 2668 2669 static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns) 2670 { 2671 Mat_IS *matis = (Mat_IS*)A->data; 2672 PetscInt nr,nl,len,i; 2673 PetscInt *lrows; 2674 2675 PetscFunctionBegin; 2676 if (PetscUnlikelyDebug(columns || diag != 0. || (x && b))) { 2677 PetscBool cong; 2678 2679 PetscCall(PetscLayoutCompare(A->rmap,A->cmap,&cong)); 2680 cong = (PetscBool)(cong && matis->sf == matis->csf); 2681 PetscCheck(cong || !columns,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Columns can be zeroed if and only if A->rmap and A->cmap are congruent and the l2g maps are the same for MATIS"); 2682 PetscCheckFalse(!cong && diag != 0.,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Nonzero diagonal value supported if and only if A->rmap and A->cmap are congruent and the l2g maps are the same for MATIS"); 2683 PetscCheckFalse(!cong && x && b,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent, and the l2g maps be the same"); 2684 } 2685 /* get locally owned rows */ 2686 PetscCall(PetscLayoutMapLocal(A->rmap,n,rows,&len,&lrows,NULL)); 2687 /* fix right hand side if needed */ 2688 if (x && b) { 2689 const PetscScalar *xx; 2690 PetscScalar *bb; 2691 2692 PetscCall(VecGetArrayRead(x, &xx)); 2693 PetscCall(VecGetArray(b, &bb)); 2694 for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 2695 PetscCall(VecRestoreArrayRead(x, &xx)); 2696 PetscCall(VecRestoreArray(b, &bb)); 2697 } 2698 /* get rows associated to the local matrices */ 2699 PetscCall(MatGetSize(matis->A,&nl,NULL)); 2700 PetscCall(PetscArrayzero(matis->sf_leafdata,nl)); 2701 PetscCall(PetscArrayzero(matis->sf_rootdata,A->rmap->n)); 2702 for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 2703 PetscCall(PetscFree(lrows)); 2704 PetscCall(PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata,MPI_REPLACE)); 2705 PetscCall(PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata,MPI_REPLACE)); 2706 PetscCall(PetscMalloc1(nl,&lrows)); 2707 for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 2708 PetscCall(MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns)); 2709 PetscCall(PetscFree(lrows)); 2710 PetscFunctionReturn(0); 2711 } 2712 2713 static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2714 { 2715 PetscFunctionBegin; 2716 PetscCall(MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE)); 2717 PetscFunctionReturn(0); 2718 } 2719 2720 static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2721 { 2722 PetscFunctionBegin; 2723 PetscCall(MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE)); 2724 PetscFunctionReturn(0); 2725 } 2726 2727 static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 2728 { 2729 Mat_IS *is = (Mat_IS*)A->data; 2730 2731 PetscFunctionBegin; 2732 PetscCall(MatAssemblyBegin(is->A,type)); 2733 PetscFunctionReturn(0); 2734 } 2735 2736 static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 2737 { 2738 Mat_IS *is = (Mat_IS*)A->data; 2739 2740 PetscFunctionBegin; 2741 PetscCall(MatAssemblyEnd(is->A,type)); 2742 /* fix for local empty rows/cols */ 2743 if (is->locempty && type == MAT_FINAL_ASSEMBLY) { 2744 Mat newlA; 2745 ISLocalToGlobalMapping rl2g,cl2g; 2746 IS nzr,nzc; 2747 PetscInt nr,nc,nnzr,nnzc; 2748 PetscBool lnewl2g,newl2g; 2749 2750 PetscCall(MatGetSize(is->A,&nr,&nc)); 2751 PetscCall(MatFindNonzeroRowsOrCols_Basic(is->A,PETSC_FALSE,PETSC_SMALL,&nzr)); 2752 if (!nzr) { 2753 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)is->A),nr,0,1,&nzr)); 2754 } 2755 PetscCall(MatFindNonzeroRowsOrCols_Basic(is->A,PETSC_TRUE,PETSC_SMALL,&nzc)); 2756 if (!nzc) { 2757 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)is->A),nc,0,1,&nzc)); 2758 } 2759 PetscCall(ISGetSize(nzr,&nnzr)); 2760 PetscCall(ISGetSize(nzc,&nnzc)); 2761 if (nnzr != nr || nnzc != nc) { /* need new global l2g map */ 2762 lnewl2g = PETSC_TRUE; 2763 PetscCallMPI(MPI_Allreduce(&lnewl2g,&newl2g,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)A))); 2764 2765 /* extract valid submatrix */ 2766 PetscCall(MatCreateSubMatrix(is->A,nzr,nzc,MAT_INITIAL_MATRIX,&newlA)); 2767 } else { /* local matrix fully populated */ 2768 lnewl2g = PETSC_FALSE; 2769 PetscCallMPI(MPI_Allreduce(&lnewl2g,&newl2g,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)A))); 2770 PetscCall(PetscObjectReference((PetscObject)is->A)); 2771 newlA = is->A; 2772 } 2773 2774 /* attach new global l2g map if needed */ 2775 if (newl2g) { 2776 IS zr,zc; 2777 const PetscInt *ridxs,*cidxs,*zridxs,*zcidxs; 2778 PetscInt *nidxs,i; 2779 2780 PetscCall(ISComplement(nzr,0,nr,&zr)); 2781 PetscCall(ISComplement(nzc,0,nc,&zc)); 2782 PetscCall(PetscMalloc1(PetscMax(nr,nc),&nidxs)); 2783 PetscCall(ISLocalToGlobalMappingGetIndices(is->rmapping,&ridxs)); 2784 PetscCall(ISLocalToGlobalMappingGetIndices(is->cmapping,&cidxs)); 2785 PetscCall(ISGetIndices(zr,&zridxs)); 2786 PetscCall(ISGetIndices(zc,&zcidxs)); 2787 PetscCall(ISGetLocalSize(zr,&nnzr)); 2788 PetscCall(ISGetLocalSize(zc,&nnzc)); 2789 2790 PetscCall(PetscArraycpy(nidxs,ridxs,nr)); 2791 for (i = 0; i < nnzr; i++) nidxs[zridxs[i]] = -1; 2792 PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,nr,nidxs,PETSC_COPY_VALUES,&rl2g)); 2793 PetscCall(PetscArraycpy(nidxs,cidxs,nc)); 2794 for (i = 0; i < nnzc; i++) nidxs[zcidxs[i]] = -1; 2795 PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,nc,nidxs,PETSC_COPY_VALUES,&cl2g)); 2796 2797 PetscCall(ISRestoreIndices(zr,&zridxs)); 2798 PetscCall(ISRestoreIndices(zc,&zcidxs)); 2799 PetscCall(ISLocalToGlobalMappingRestoreIndices(is->rmapping,&ridxs)); 2800 PetscCall(ISLocalToGlobalMappingRestoreIndices(is->cmapping,&cidxs)); 2801 PetscCall(ISDestroy(&nzr)); 2802 PetscCall(ISDestroy(&nzc)); 2803 PetscCall(ISDestroy(&zr)); 2804 PetscCall(ISDestroy(&zc)); 2805 PetscCall(PetscFree(nidxs)); 2806 PetscCall(MatSetLocalToGlobalMapping(A,rl2g,cl2g)); 2807 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 2808 PetscCall(ISLocalToGlobalMappingDestroy(&cl2g)); 2809 } 2810 PetscCall(MatISSetLocalMat(A,newlA)); 2811 PetscCall(MatDestroy(&newlA)); 2812 PetscCall(ISDestroy(&nzr)); 2813 PetscCall(ISDestroy(&nzc)); 2814 is->locempty = PETSC_FALSE; 2815 } 2816 PetscFunctionReturn(0); 2817 } 2818 2819 static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 2820 { 2821 Mat_IS *is = (Mat_IS*)mat->data; 2822 2823 PetscFunctionBegin; 2824 *local = is->A; 2825 PetscFunctionReturn(0); 2826 } 2827 2828 static PetscErrorCode MatISRestoreLocalMat_IS(Mat mat,Mat *local) 2829 { 2830 PetscFunctionBegin; 2831 *local = NULL; 2832 PetscFunctionReturn(0); 2833 } 2834 2835 /*@ 2836 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 2837 2838 Input Parameter: 2839 . mat - the matrix 2840 2841 Output Parameter: 2842 . local - the local matrix 2843 2844 Level: advanced 2845 2846 Notes: 2847 This can be called if you have precomputed the nonzero structure of the 2848 matrix and want to provide it to the inner matrix object to improve the performance 2849 of the MatSetValues() operation. 2850 2851 Call MatISRestoreLocalMat() when finished with the local matrix. 2852 2853 .seealso: MATIS 2854 @*/ 2855 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 2856 { 2857 PetscFunctionBegin; 2858 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2859 PetscValidPointer(local,2); 2860 PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local)); 2861 PetscFunctionReturn(0); 2862 } 2863 2864 /*@ 2865 MatISRestoreLocalMat - Restores the local matrix obtained with MatISGetLocalMat() 2866 2867 Input Parameter: 2868 . mat - the matrix 2869 2870 Output Parameter: 2871 . local - the local matrix 2872 2873 Level: advanced 2874 2875 .seealso: MATIS 2876 @*/ 2877 PetscErrorCode MatISRestoreLocalMat(Mat mat,Mat *local) 2878 { 2879 PetscFunctionBegin; 2880 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2881 PetscValidPointer(local,2); 2882 PetscUseMethod(mat,"MatISRestoreLocalMat_C",(Mat,Mat*),(mat,local)); 2883 PetscFunctionReturn(0); 2884 } 2885 2886 static PetscErrorCode MatISSetLocalMatType_IS(Mat mat,MatType mtype) 2887 { 2888 Mat_IS *is = (Mat_IS*)mat->data; 2889 2890 PetscFunctionBegin; 2891 if (is->A) { 2892 PetscCall(MatSetType(is->A,mtype)); 2893 } 2894 PetscCall(PetscFree(is->lmattype)); 2895 PetscCall(PetscStrallocpy(mtype,&is->lmattype)); 2896 PetscFunctionReturn(0); 2897 } 2898 2899 /*@ 2900 MatISSetLocalMatType - Specifies the type of local matrix 2901 2902 Input Parameters: 2903 + mat - the matrix 2904 - mtype - the local matrix type 2905 2906 Output Parameter: 2907 2908 Level: advanced 2909 2910 .seealso: MATIS, MatSetType(), MatType 2911 @*/ 2912 PetscErrorCode MatISSetLocalMatType(Mat mat,MatType mtype) 2913 { 2914 PetscFunctionBegin; 2915 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2916 PetscUseMethod(mat,"MatISSetLocalMatType_C",(Mat,MatType),(mat,mtype)); 2917 PetscFunctionReturn(0); 2918 } 2919 2920 static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 2921 { 2922 Mat_IS *is = (Mat_IS*)mat->data; 2923 PetscInt nrows,ncols,orows,ocols; 2924 MatType mtype,otype; 2925 PetscBool sametype = PETSC_TRUE; 2926 2927 PetscFunctionBegin; 2928 if (is->A && !is->islocalref) { 2929 PetscCall(MatGetSize(is->A,&orows,&ocols)); 2930 PetscCall(MatGetSize(local,&nrows,&ncols)); 2931 PetscCheckFalse(orows != nrows || ocols != ncols,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %" PetscInt_FMT "x%" PetscInt_FMT " (you passed a %" PetscInt_FMT "x%" PetscInt_FMT " matrix)",orows,ocols,nrows,ncols); 2932 PetscCall(MatGetType(local,&mtype)); 2933 PetscCall(MatGetType(is->A,&otype)); 2934 PetscCall(PetscStrcmp(mtype,otype,&sametype)); 2935 } 2936 PetscCall(PetscObjectReference((PetscObject)local)); 2937 PetscCall(MatDestroy(&is->A)); 2938 is->A = local; 2939 PetscCall(MatGetType(is->A,&mtype)); 2940 PetscCall(MatISSetLocalMatType(mat,mtype)); 2941 if (!sametype && !is->islocalref) { 2942 PetscCall(MatISSetUpScatters_Private(mat)); 2943 } 2944 PetscFunctionReturn(0); 2945 } 2946 2947 /*@ 2948 MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 2949 2950 Collective on Mat 2951 2952 Input Parameters: 2953 + mat - the matrix 2954 - local - the local matrix 2955 2956 Output Parameter: 2957 2958 Level: advanced 2959 2960 Notes: 2961 This can be called if you have precomputed the local matrix and 2962 want to provide it to the matrix object MATIS. 2963 2964 .seealso: MATIS 2965 @*/ 2966 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 2967 { 2968 PetscFunctionBegin; 2969 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2970 PetscValidHeaderSpecific(local,MAT_CLASSID,2); 2971 PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local)); 2972 PetscFunctionReturn(0); 2973 } 2974 2975 static PetscErrorCode MatZeroEntries_IS(Mat A) 2976 { 2977 Mat_IS *a = (Mat_IS*)A->data; 2978 2979 PetscFunctionBegin; 2980 PetscCall(MatZeroEntries(a->A)); 2981 PetscFunctionReturn(0); 2982 } 2983 2984 static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 2985 { 2986 Mat_IS *is = (Mat_IS*)A->data; 2987 2988 PetscFunctionBegin; 2989 PetscCall(MatScale(is->A,a)); 2990 PetscFunctionReturn(0); 2991 } 2992 2993 static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 2994 { 2995 Mat_IS *is = (Mat_IS*)A->data; 2996 2997 PetscFunctionBegin; 2998 /* get diagonal of the local matrix */ 2999 PetscCall(MatGetDiagonal(is->A,is->y)); 3000 3001 /* scatter diagonal back into global vector */ 3002 PetscCall(VecSet(v,0)); 3003 PetscCall(VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE)); 3004 PetscCall(VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE)); 3005 PetscFunctionReturn(0); 3006 } 3007 3008 static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 3009 { 3010 Mat_IS *a = (Mat_IS*)A->data; 3011 3012 PetscFunctionBegin; 3013 PetscCall(MatSetOption(a->A,op,flg)); 3014 PetscFunctionReturn(0); 3015 } 3016 3017 static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 3018 { 3019 Mat_IS *y = (Mat_IS*)Y->data; 3020 Mat_IS *x; 3021 3022 PetscFunctionBegin; 3023 if (PetscDefined(USE_DEBUG)) { 3024 PetscBool ismatis; 3025 PetscCall(PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis)); 3026 PetscCheck(ismatis,PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 3027 } 3028 x = (Mat_IS*)X->data; 3029 PetscCall(MatAXPY(y->A,a,x->A,str)); 3030 PetscFunctionReturn(0); 3031 } 3032 3033 static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 3034 { 3035 Mat lA; 3036 Mat_IS *matis = (Mat_IS*)(A->data); 3037 ISLocalToGlobalMapping rl2g,cl2g; 3038 IS is; 3039 const PetscInt *rg,*rl; 3040 PetscInt nrg; 3041 PetscInt N,M,nrl,i,*idxs; 3042 3043 PetscFunctionBegin; 3044 PetscCall(ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg)); 3045 PetscCall(ISGetLocalSize(row,&nrl)); 3046 PetscCall(ISGetIndices(row,&rl)); 3047 PetscCall(ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg)); 3048 if (PetscDefined(USE_DEBUG)) { 3049 for (i=0; i<nrl; i++) PetscCheck(rl[i]<nrg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row index %" PetscInt_FMT " -> %" PetscInt_FMT " greater then maximum possible %" PetscInt_FMT,i,rl[i],nrg); 3050 } 3051 PetscCall(PetscMalloc1(nrg,&idxs)); 3052 /* map from [0,nrl) to row */ 3053 for (i=0;i<nrl;i++) idxs[i] = rl[i]; 3054 for (i=nrl;i<nrg;i++) idxs[i] = -1; 3055 PetscCall(ISRestoreIndices(row,&rl)); 3056 PetscCall(ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg)); 3057 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is)); 3058 PetscCall(ISLocalToGlobalMappingCreateIS(is,&rl2g)); 3059 PetscCall(ISDestroy(&is)); 3060 /* compute new l2g map for columns */ 3061 if (col != row || matis->rmapping != matis->cmapping || matis->A->rmap->mapping != matis->A->cmap->mapping) { 3062 const PetscInt *cg,*cl; 3063 PetscInt ncg; 3064 PetscInt ncl; 3065 3066 PetscCall(ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg)); 3067 PetscCall(ISGetLocalSize(col,&ncl)); 3068 PetscCall(ISGetIndices(col,&cl)); 3069 PetscCall(ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg)); 3070 if (PetscDefined(USE_DEBUG)) { 3071 for (i=0; i<ncl; i++) PetscCheck(cl[i]<ncg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local column index %" PetscInt_FMT " -> %" PetscInt_FMT " greater then maximum possible %" PetscInt_FMT,i,cl[i],ncg); 3072 } 3073 PetscCall(PetscMalloc1(ncg,&idxs)); 3074 /* map from [0,ncl) to col */ 3075 for (i=0;i<ncl;i++) idxs[i] = cl[i]; 3076 for (i=ncl;i<ncg;i++) idxs[i] = -1; 3077 PetscCall(ISRestoreIndices(col,&cl)); 3078 PetscCall(ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg)); 3079 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is)); 3080 PetscCall(ISLocalToGlobalMappingCreateIS(is,&cl2g)); 3081 PetscCall(ISDestroy(&is)); 3082 } else { 3083 PetscCall(PetscObjectReference((PetscObject)rl2g)); 3084 cl2g = rl2g; 3085 } 3086 /* create the MATIS submatrix */ 3087 PetscCall(MatGetSize(A,&M,&N)); 3088 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),submat)); 3089 PetscCall(MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N)); 3090 PetscCall(MatSetType(*submat,MATIS)); 3091 matis = (Mat_IS*)((*submat)->data); 3092 matis->islocalref = PETSC_TRUE; 3093 PetscCall(MatSetLocalToGlobalMapping(*submat,rl2g,cl2g)); 3094 PetscCall(MatISGetLocalMat(A,&lA)); 3095 PetscCall(MatISSetLocalMat(*submat,lA)); 3096 PetscCall(MatSetUp(*submat)); 3097 PetscCall(MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY)); 3098 PetscCall(MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY)); 3099 PetscCall(ISLocalToGlobalMappingDestroy(&rl2g)); 3100 PetscCall(ISLocalToGlobalMappingDestroy(&cl2g)); 3101 3102 /* remove unsupported ops */ 3103 PetscCall(PetscMemzero((*submat)->ops,sizeof(struct _MatOps))); 3104 (*submat)->ops->destroy = MatDestroy_IS; 3105 (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 3106 (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 3107 (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 3108 (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 3109 PetscFunctionReturn(0); 3110 } 3111 3112 static PetscErrorCode MatSetFromOptions_IS(PetscOptionItems *PetscOptionsObject, Mat A) 3113 { 3114 Mat_IS *a = (Mat_IS*)A->data; 3115 char type[256]; 3116 PetscBool flg; 3117 3118 PetscFunctionBegin; 3119 PetscOptionsHeadBegin(PetscOptionsObject,"MATIS options"); 3120 PetscCall(PetscOptionsBool("-matis_fixempty","Fix local matrices in case of empty local rows/columns","MatISFixLocalEmpty",a->locempty,&a->locempty,NULL)); 3121 PetscCall(PetscOptionsBool("-matis_storel2l","Store local-to-local matrices generated from PtAP operations","MatISStoreL2L",a->storel2l,&a->storel2l,NULL)); 3122 PetscCall(PetscOptionsFList("-matis_localmat_type","Matrix type","MatISSetLocalMatType",MatList,a->lmattype,type,256,&flg)); 3123 if (flg) { 3124 PetscCall(MatISSetLocalMatType(A,type)); 3125 } 3126 if (a->A) { 3127 PetscCall(MatSetFromOptions(a->A)); 3128 } 3129 PetscOptionsHeadEnd(); 3130 PetscFunctionReturn(0); 3131 } 3132 3133 /*@ 3134 MatCreateIS - Creates a "process" unassembled matrix, assembled on each 3135 process but not across processes. 3136 3137 Input Parameters: 3138 + comm - MPI communicator that will share the matrix 3139 . bs - block size of the matrix 3140 . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 3141 . rmap - local to global map for rows 3142 - cmap - local to global map for cols 3143 3144 Output Parameter: 3145 . A - the resulting matrix 3146 3147 Level: advanced 3148 3149 Notes: 3150 See MATIS for more details. 3151 m and n are NOT related to the size of the map; they represent the size of the local parts of the distributed vectors 3152 used in MatMult operations. The sizes of rmap and cmap define the size of the local matrices. 3153 If rmap (cmap) is NULL, then the local row (column) spaces matches the global space. 3154 3155 .seealso: MATIS, MatSetLocalToGlobalMapping() 3156 @*/ 3157 PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 3158 { 3159 PetscFunctionBegin; 3160 PetscCall(MatCreate(comm,A)); 3161 PetscCall(MatSetSizes(*A,m,n,M,N)); 3162 if (bs > 0) { 3163 PetscCall(MatSetBlockSize(*A,bs)); 3164 } 3165 PetscCall(MatSetType(*A,MATIS)); 3166 PetscCall(MatSetLocalToGlobalMapping(*A,rmap,cmap)); 3167 PetscFunctionReturn(0); 3168 } 3169 3170 static PetscErrorCode MatHasOperation_IS(Mat A, MatOperation op, PetscBool *has) 3171 { 3172 Mat_IS *a = (Mat_IS*)A->data; 3173 3174 PetscFunctionBegin; 3175 *has = PETSC_FALSE; 3176 if (!((void**)A->ops)[op]) PetscFunctionReturn(0); 3177 PetscCall(MatHasOperation(a->A,op,has)); 3178 PetscFunctionReturn(0); 3179 } 3180 3181 static PetscErrorCode MatSetValuesCOO_IS(Mat A,const PetscScalar v[],InsertMode imode) 3182 { 3183 Mat_IS *a = (Mat_IS*)A->data; 3184 3185 PetscFunctionBegin; 3186 PetscCall(MatSetValuesCOO(a->A,v,imode)); 3187 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 3188 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 3189 PetscFunctionReturn(0); 3190 } 3191 3192 static PetscErrorCode MatSetPreallocationCOOLocal_IS(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[]) 3193 { 3194 Mat_IS *a = (Mat_IS*)A->data; 3195 3196 PetscFunctionBegin; 3197 PetscCheck(a->A,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to provide l2g map first via MatSetLocalToGlobalMapping"); 3198 if (a->A->rmap->mapping || a->A->cmap->mapping) { 3199 PetscCall(MatSetPreallocationCOOLocal(a->A,ncoo,coo_i,coo_j)); 3200 } else { 3201 PetscCall(MatSetPreallocationCOO(a->A,ncoo,coo_i,coo_j)); 3202 } 3203 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetValuesCOO_C",MatSetValuesCOO_IS)); 3204 A->preallocated = PETSC_TRUE; 3205 PetscFunctionReturn(0); 3206 } 3207 3208 static PetscErrorCode MatSetPreallocationCOO_IS(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 3209 { 3210 Mat_IS *a = (Mat_IS*)A->data; 3211 PetscInt *coo_il, *coo_jl, incoo; 3212 3213 PetscFunctionBegin; 3214 PetscCheck(a->A,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to provide l2g map first via MatSetLocalToGlobalMapping"); 3215 PetscCheck(ncoo <= PETSC_MAX_INT,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support",ncoo); 3216 PetscCall(PetscMalloc2(ncoo,&coo_il,ncoo,&coo_jl)); 3217 PetscCall(ISGlobalToLocalMappingApply(a->rmapping,IS_GTOLM_MASK,ncoo,coo_i,&incoo,coo_il)); 3218 PetscCall(ISGlobalToLocalMappingApply(a->cmapping,IS_GTOLM_MASK,ncoo,coo_j,&incoo,coo_jl)); 3219 PetscCall(MatSetPreallocationCOO(a->A,ncoo,coo_il,coo_jl)); 3220 PetscCall(PetscFree2(coo_il,coo_jl)); 3221 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetValuesCOO_C",MatSetValuesCOO_IS)); 3222 A->preallocated = PETSC_TRUE; 3223 PetscFunctionReturn(0); 3224 } 3225 3226 /*@ 3227 MatISGetLocalToGlobalMapping - Gets the local-to-global numbering of the MATIS object 3228 3229 Not Collective 3230 3231 Input Parameter: 3232 . A - the matrix 3233 3234 Output Parameters: 3235 + rmapping - row mapping 3236 - cmapping - column mapping 3237 3238 Notes: The returned map can be different from the one used to construct the MATIS object, since it will not contain negative or repeated indices. 3239 3240 Level: advanced 3241 3242 .seealso: MatSetLocalToGlobalMapping() 3243 @*/ 3244 PetscErrorCode MatISGetLocalToGlobalMapping(Mat A,ISLocalToGlobalMapping *rmapping,ISLocalToGlobalMapping *cmapping) 3245 { 3246 PetscFunctionBegin; 3247 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3248 PetscValidType(A,1); 3249 if (rmapping) PetscValidPointer(rmapping,2); 3250 if (cmapping) PetscValidPointer(cmapping,3); 3251 PetscUseMethod(A,"MatISGetLocalToGlobalMapping_C",(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*),(A,rmapping,cmapping)); 3252 PetscFunctionReturn(0); 3253 } 3254 3255 static PetscErrorCode MatISGetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping *r, ISLocalToGlobalMapping *c) 3256 { 3257 Mat_IS *a = (Mat_IS*)A->data; 3258 3259 PetscFunctionBegin; 3260 if (r) *r = a->rmapping; 3261 if (c) *c = a->cmapping; 3262 PetscFunctionReturn(0); 3263 } 3264 3265 /*MC 3266 MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 3267 This stores the matrices in globally unassembled form. Each processor assembles only its local Neumann problem and the parallel matrix vector 3268 product is handled "implicitly". 3269 3270 Options Database Keys: 3271 + -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 3272 . -matis_fixempty - Fixes local matrices in case of empty local rows/columns. 3273 - -matis_storel2l - stores the local-to-local operators generated by the Galerkin process of MatPtAP(). 3274 3275 Notes: 3276 Options prefix for the inner matrix are given by -is_mat_xxx 3277 3278 You must call MatSetLocalToGlobalMapping() before using this matrix type. 3279 3280 You can do matrix preallocation on the local matrix after you obtain it with 3281 MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 3282 3283 Level: advanced 3284 3285 .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 3286 3287 M*/ 3288 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 3289 { 3290 Mat_IS *a; 3291 3292 PetscFunctionBegin; 3293 PetscCall(PetscNewLog(A,&a)); 3294 PetscCall(PetscStrallocpy(MATAIJ,&a->lmattype)); 3295 A->data = (void*)a; 3296 3297 /* matrix ops */ 3298 PetscCall(PetscMemzero(A->ops,sizeof(struct _MatOps))); 3299 A->ops->mult = MatMult_IS; 3300 A->ops->multadd = MatMultAdd_IS; 3301 A->ops->multtranspose = MatMultTranspose_IS; 3302 A->ops->multtransposeadd = MatMultTransposeAdd_IS; 3303 A->ops->destroy = MatDestroy_IS; 3304 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 3305 A->ops->setvalues = MatSetValues_IS; 3306 A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 3307 A->ops->setvalueslocal = MatSetValuesLocal_IS; 3308 A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 3309 A->ops->zerorows = MatZeroRows_IS; 3310 A->ops->zerorowscolumns = MatZeroRowsColumns_IS; 3311 A->ops->assemblybegin = MatAssemblyBegin_IS; 3312 A->ops->assemblyend = MatAssemblyEnd_IS; 3313 A->ops->view = MatView_IS; 3314 A->ops->zeroentries = MatZeroEntries_IS; 3315 A->ops->scale = MatScale_IS; 3316 A->ops->getdiagonal = MatGetDiagonal_IS; 3317 A->ops->setoption = MatSetOption_IS; 3318 A->ops->ishermitian = MatIsHermitian_IS; 3319 A->ops->issymmetric = MatIsSymmetric_IS; 3320 A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS; 3321 A->ops->duplicate = MatDuplicate_IS; 3322 A->ops->missingdiagonal = MatMissingDiagonal_IS; 3323 A->ops->copy = MatCopy_IS; 3324 A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 3325 A->ops->createsubmatrix = MatCreateSubMatrix_IS; 3326 A->ops->axpy = MatAXPY_IS; 3327 A->ops->diagonalset = MatDiagonalSet_IS; 3328 A->ops->shift = MatShift_IS; 3329 A->ops->transpose = MatTranspose_IS; 3330 A->ops->getinfo = MatGetInfo_IS; 3331 A->ops->diagonalscale = MatDiagonalScale_IS; 3332 A->ops->setfromoptions = MatSetFromOptions_IS; 3333 A->ops->setup = MatSetUp_IS; 3334 A->ops->hasoperation = MatHasOperation_IS; 3335 3336 /* special MATIS functions */ 3337 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMatType_C",MatISSetLocalMatType_IS)); 3338 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS)); 3339 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISRestoreLocalMat_C",MatISRestoreLocalMat_IS)); 3340 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS)); 3341 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatConvert_IS_XAIJ)); 3342 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS)); 3343 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISStoreL2L_C",MatISStoreL2L_IS)); 3344 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISFixLocalEmpty_C",MatISFixLocalEmpty_IS)); 3345 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalToGlobalMapping_C",MatISGetLocalToGlobalMapping_IS)); 3346 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_mpiaij_C",MatConvert_IS_XAIJ)); 3347 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_mpibaij_C",MatConvert_IS_XAIJ)); 3348 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_mpisbaij_C",MatConvert_IS_XAIJ)); 3349 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_seqaij_C",MatConvert_IS_XAIJ)); 3350 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_seqbaij_C",MatConvert_IS_XAIJ)); 3351 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_seqsbaij_C",MatConvert_IS_XAIJ)); 3352 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_is_aij_C",MatConvert_IS_XAIJ)); 3353 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",MatSetPreallocationCOOLocal_IS)); 3354 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSetPreallocationCOO_C",MatSetPreallocationCOO_IS)); 3355 PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATIS)); 3356 PetscFunctionReturn(0); 3357 } 3358