1 2 #include "matnestimpl.h" /*I "petscmat.h" I*/ 3 4 static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]); 5 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left); 6 7 /* private functions */ 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatNestGetSizes_Private" 10 static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N) 11 { 12 Mat_Nest *bA = (Mat_Nest*)A->data; 13 PetscInt i,j; 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 *m = *n = *M = *N = 0; 18 for (i=0; i<bA->nr; i++) { /* rows */ 19 PetscInt sm,sM; 20 ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr); 21 ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr); 22 *m += sm; 23 *M += sM; 24 } 25 for (j=0; j<bA->nc; j++) { /* cols */ 26 PetscInt sn,sN; 27 ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr); 28 ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr); 29 *n += sn; 30 *N += sN; 31 } 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2" 37 PETSC_UNUSED 38 static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y) 39 { 40 PetscBool isANest, isxNest, isyNest; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 isANest = isxNest = isyNest = PETSC_FALSE; 45 ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr); 46 ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr); 47 ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr); 48 49 if (isANest && isxNest && isyNest) { 50 PetscFunctionReturn(0); 51 } else { 52 SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations"); 53 PetscFunctionReturn(0); 54 } 55 PetscFunctionReturn(0); 56 } 57 58 /* 59 This function is called every time we insert a sub matrix the Nest. 60 It ensures that every Nest along row r, and col c has the same dimensions 61 as the submat being inserted. 62 */ 63 #undef __FUNCT__ 64 #define __FUNCT__ "PETSc_MatNest_CheckConsistency" 65 PETSC_UNUSED 66 static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c) 67 { 68 Mat_Nest *b = (Mat_Nest*)A->data; 69 PetscInt i,j; 70 PetscInt nr,nc; 71 PetscInt sM,sN,M,N; 72 Mat mat; 73 PetscErrorCode ierr; 74 75 PetscFunctionBegin; 76 nr = b->nr; 77 nc = b->nc; 78 ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr); 79 for (i=0; i<nr; i++) { 80 mat = b->m[i][c]; 81 if (mat) { 82 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 83 /* Check columns match */ 84 if (sN != N) { 85 SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N ); 86 } 87 } 88 } 89 90 for (j=0; j<nc; j++) { 91 mat = b->m[r][j]; 92 if (mat) { 93 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 94 /* Check rows match */ 95 if (sM != M) { 96 SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M ); 97 } 98 } 99 } 100 PetscFunctionReturn(0); 101 } 102 103 /* 104 Checks if the row/col's contain a complete set of IS's. 105 Once they are filled, the offsets are computed. This saves having to update 106 the index set entries each time we insert something new. 107 */ 108 #undef __FUNCT__ 109 #define __FUNCT__ "PETSc_MatNest_UpdateStructure" 110 PETSC_UNUSED 111 static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx) 112 { 113 Mat_Nest *b = (Mat_Nest*)A->data; 114 PetscInt m,n,i,j; 115 PetscBool fullrow,fullcol; 116 PetscErrorCode ierr; 117 118 PetscFunctionBegin; 119 ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr); 120 b->row_len[ridx] = m; 121 b->col_len[cidx] = n; 122 123 fullrow = PETSC_TRUE; 124 for (i=0; i<b->nr; i++) { 125 if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; } 126 } 127 if ( (fullrow) && (!b->isglobal.row[0]) ){ 128 PetscInt cnt = 0; 129 for (i=0; i<b->nr; i++) { 130 ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr); 131 cnt = cnt + b->row_len[i]; 132 } 133 } 134 135 fullcol = PETSC_TRUE; 136 for (j=0; j<b->nc; j++) { 137 if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; } 138 } 139 if( (fullcol) && (!b->isglobal.col[0]) ){ 140 PetscInt cnt = 0; 141 for (j=0; j<b->nc; j++) { 142 ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr); 143 cnt = cnt + b->col_len[j]; 144 } 145 } 146 PetscFunctionReturn(0); 147 } 148 149 /* operations */ 150 #undef __FUNCT__ 151 #define __FUNCT__ "MatMult_Nest" 152 static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 153 { 154 Mat_Nest *bA = (Mat_Nest*)A->data; 155 Vec *bx = bA->right,*by = bA->left; 156 PetscInt i,j,nr = bA->nr,nc = bA->nc; 157 PetscErrorCode ierr; 158 159 PetscFunctionBegin; 160 for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 161 for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 162 for (i=0; i<nr; i++) { 163 ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 164 for (j=0; j<nc; j++) { 165 if (!bA->m[i][j]) continue; 166 /* y[i] <- y[i] + A[i][j] * x[j] */ 167 ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 168 } 169 } 170 for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 171 for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 172 PetscFunctionReturn(0); 173 } 174 175 #undef __FUNCT__ 176 #define __FUNCT__ "MatMultTranspose_Nest" 177 static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 178 { 179 Mat_Nest *bA = (Mat_Nest*)A->data; 180 Vec *bx = bA->left,*by = bA->right; 181 PetscInt i,j,nr = bA->nr,nc = bA->nc; 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 186 for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 187 for (j=0; j<nc; j++) { 188 ierr = VecZeroEntries(by[j]);CHKERRQ(ierr); 189 for (i=0; i<nr; i++) { 190 if (!bA->m[j][i]) continue; 191 /* y[j] <- y[j] + (A[i][j])^T * x[i] */ 192 ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr); 193 } 194 } 195 for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 196 for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 197 PetscFunctionReturn(0); 198 } 199 200 #undef __FUNCT__ 201 #define __FUNCT__ "MatNestDestroyISList" 202 static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 203 { 204 PetscErrorCode ierr; 205 IS *lst = *list; 206 PetscInt i; 207 208 PetscFunctionBegin; 209 if (!lst) PetscFunctionReturn(0); 210 for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);} 211 ierr = PetscFree(lst);CHKERRQ(ierr); 212 *list = PETSC_NULL; 213 PetscFunctionReturn(0); 214 } 215 216 #undef __FUNCT__ 217 #define __FUNCT__ "MatDestroy_Nest" 218 static PetscErrorCode MatDestroy_Nest(Mat A) 219 { 220 Mat_Nest *vs = (Mat_Nest*)A->data; 221 PetscInt i,j; 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 /* release the matrices and the place holders */ 226 ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 227 ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 228 ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 229 ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 230 231 ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 232 ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 233 234 ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 235 236 /* release the matrices and the place holders */ 237 if (vs->m) { 238 for (i=0; i<vs->nr; i++) { 239 for (j=0; j<vs->nc; j++) { 240 ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr); 241 } 242 ierr = PetscFree( vs->m[i] );CHKERRQ(ierr); 243 } 244 ierr = PetscFree(vs->m);CHKERRQ(ierr); 245 } 246 ierr = PetscFree(A->data);CHKERRQ(ierr); 247 248 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "",0);CHKERRQ(ierr); 249 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "",0);CHKERRQ(ierr); 250 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "",0);CHKERRQ(ierr); 251 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "",0);CHKERRQ(ierr); 252 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "",0);CHKERRQ(ierr); 253 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "",0);CHKERRQ(ierr); 254 PetscFunctionReturn(0); 255 } 256 257 #undef __FUNCT__ 258 #define __FUNCT__ "MatAssemblyBegin_Nest" 259 static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 260 { 261 Mat_Nest *vs = (Mat_Nest*)A->data; 262 PetscInt i,j; 263 PetscErrorCode ierr; 264 265 PetscFunctionBegin; 266 for (i=0; i<vs->nr; i++) { 267 for (j=0; j<vs->nc; j++) { 268 if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); } 269 } 270 } 271 PetscFunctionReturn(0); 272 } 273 274 #undef __FUNCT__ 275 #define __FUNCT__ "MatAssemblyEnd_Nest" 276 static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 277 { 278 Mat_Nest *vs = (Mat_Nest*)A->data; 279 PetscInt i,j; 280 PetscErrorCode ierr; 281 282 PetscFunctionBegin; 283 for (i=0; i<vs->nr; i++) { 284 for (j=0; j<vs->nc; j++) { 285 if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); } 286 } 287 } 288 PetscFunctionReturn(0); 289 } 290 291 #undef __FUNCT__ 292 #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 293 static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 294 { 295 PetscErrorCode ierr; 296 Mat_Nest *vs = (Mat_Nest*)A->data; 297 PetscInt j; 298 Mat sub; 299 300 PetscFunctionBegin; 301 sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */ 302 for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 303 if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 304 *B = sub; 305 PetscFunctionReturn(0); 306 } 307 308 #undef __FUNCT__ 309 #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 310 static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 311 { 312 PetscErrorCode ierr; 313 Mat_Nest *vs = (Mat_Nest*)A->data; 314 PetscInt i; 315 Mat sub; 316 317 PetscFunctionBegin; 318 sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */ 319 for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 320 if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 321 *B = sub; 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "MatNestFindIS" 327 static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 328 { 329 PetscErrorCode ierr; 330 PetscInt i; 331 PetscBool flg; 332 333 PetscFunctionBegin; 334 PetscValidPointer(list,3); 335 PetscValidHeaderSpecific(is,IS_CLASSID,4); 336 PetscValidIntPointer(found,5); 337 *found = -1; 338 for (i=0; i<n; i++) { 339 if (!list[i]) continue; 340 ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 341 if (flg) { 342 *found = i; 343 PetscFunctionReturn(0); 344 } 345 } 346 SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set"); 347 PetscFunctionReturn(0); 348 } 349 350 #undef __FUNCT__ 351 #define __FUNCT__ "MatNestGetRow" 352 /* Get a block row as a new MatNest */ 353 static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 354 { 355 Mat_Nest *vs = (Mat_Nest*)A->data; 356 Mat C; 357 char keyname[256]; 358 PetscErrorCode ierr; 359 360 PetscFunctionBegin; 361 *B = PETSC_NULL; 362 ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr); 363 ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 364 if (*B) PetscFunctionReturn(0); 365 366 ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 367 (*B)->assembled = A->assembled; 368 ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 369 ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 370 PetscFunctionReturn(0); 371 } 372 373 #undef __FUNCT__ 374 #define __FUNCT__ "MatNestFindSubMat" 375 static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 376 { 377 Mat_Nest *vs = (Mat_Nest*)A->data; 378 PetscErrorCode ierr; 379 PetscInt row,col,i; 380 IS *iscopy; 381 Mat *matcopy; 382 PetscBool same,isFullCol; 383 384 PetscFunctionBegin; 385 /* Check if full column space. This is a hack */ 386 isFullCol = PETSC_FALSE; 387 ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 388 if (same) { 389 PetscInt n,first,step; 390 ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 391 ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 392 if (A->cmap->n == n && A->cmap->rstart == first && step == 1) isFullCol = PETSC_TRUE; 393 } 394 395 if (isFullCol) { 396 PetscInt row; 397 ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 398 ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 399 } else { 400 ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 401 ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 402 *B = vs->m[row][col]; 403 } 404 PetscFunctionReturn(0); 405 } 406 407 #undef __FUNCT__ 408 #define __FUNCT__ "MatGetSubMatrix_Nest" 409 static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 410 { 411 PetscErrorCode ierr; 412 Mat_Nest *vs = (Mat_Nest*)A->data; 413 Mat sub; 414 415 PetscFunctionBegin; 416 ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 417 switch (reuse) { 418 case MAT_INITIAL_MATRIX: 419 if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 420 *B = sub; 421 break; 422 case MAT_REUSE_MATRIX: 423 if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 424 break; 425 case MAT_IGNORE_MATRIX: /* Nothing to do */ 426 break; 427 } 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 433 PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 434 { 435 PetscErrorCode ierr; 436 Mat_Nest *vs = (Mat_Nest*)A->data; 437 Mat sub; 438 439 PetscFunctionBegin; 440 ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 441 /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 442 if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 443 *B = sub; 444 PetscFunctionReturn(0); 445 } 446 447 #undef __FUNCT__ 448 #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 449 static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 450 { 451 PetscErrorCode ierr; 452 Mat_Nest *vs = (Mat_Nest*)A->data; 453 Mat sub; 454 455 PetscFunctionBegin; 456 ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 457 if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 458 if (sub) { 459 if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 460 ierr = MatDestroy(B);CHKERRQ(ierr); 461 } 462 PetscFunctionReturn(0); 463 } 464 465 #undef __FUNCT__ 466 #define __FUNCT__ "MatGetDiagonal_Nest" 467 static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 468 { 469 Mat_Nest *bA = (Mat_Nest*)A->data; 470 Vec *bdiag; 471 PetscInt i; 472 PetscErrorCode ierr; 473 474 PetscFunctionBegin; 475 /* ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); */ 476 /* ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); */ 477 ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr); 478 for (i=0; i<bA->nr; i++) { 479 if (bA->m[i][i]) { 480 ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr); 481 } else { 482 ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr); 483 } 484 } 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "MatGetDiagonal_Nest2" 490 static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v) 491 { 492 Mat_Nest *bA = (Mat_Nest*)A->data; 493 Vec diag,*bdiag; 494 VecScatter *vscat; 495 PetscInt i; 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); 500 ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr); 501 502 /* scatter diag into v */ 503 ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); 504 ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr); 505 for (i=0; i<bA->nr; i++) { 506 ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr); 507 ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 508 } 509 for (i=0; i<bA->nr; i++) { 510 ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 511 } 512 for (i=0; i<bA->nr; i++) { 513 ierr = VecScatterDestroy(&vscat[i]);CHKERRQ(ierr); 514 } 515 ierr = PetscFree(vscat);CHKERRQ(ierr); 516 ierr = VecDestroy(&diag);CHKERRQ(ierr); 517 PetscFunctionReturn(0); 518 } 519 520 #undef __FUNCT__ 521 #define __FUNCT__ "MatDiagonalScale_Nest" 522 static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 523 { 524 Mat_Nest *bA = (Mat_Nest*)A->data; 525 Vec *bl,*br; 526 PetscInt i,j; 527 PetscErrorCode ierr; 528 529 PetscFunctionBegin; 530 ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr); 531 ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr); 532 for (i=0; i<bA->nr; i++) { 533 for (j=0; j<bA->nc; j++) { 534 if (bA->m[i][j]) { 535 ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr); 536 } 537 } 538 } 539 PetscFunctionReturn(0); 540 } 541 542 #undef __FUNCT__ 543 #define __FUNCT__ "MatDiagonalScale_Nest2" 544 static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r) 545 { 546 Mat_Nest *bA = (Mat_Nest*)A->data; 547 Vec bl,br,*ble,*bre; 548 VecScatter *vscatl,*vscatr; 549 PetscInt i; 550 PetscErrorCode ierr; 551 552 PetscFunctionBegin; 553 /* scatter l,r into bl,br */ 554 ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr); 555 ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr); 556 ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr); 557 558 /* row */ 559 ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr); 560 for (i=0; i<bA->nr; i++) { 561 ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr); 562 ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 563 } 564 for (i=0; i<bA->nr; i++) { 565 ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 566 } 567 /* col */ 568 ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr); 569 for (i=0; i<bA->nc; i++) { 570 ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr); 571 ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 572 } 573 for (i=0; i<bA->nc; i++) { 574 ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 575 } 576 577 ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr); 578 579 for (i=0; i<bA->nr; i++) { 580 ierr = VecScatterDestroy(&vscatl[i]);CHKERRQ(ierr); 581 } 582 for (i=0; i<bA->nc; i++) { 583 ierr = VecScatterDestroy(&vscatr[i]);CHKERRQ(ierr); 584 } 585 ierr = PetscFree(vscatl);CHKERRQ(ierr); 586 ierr = PetscFree(vscatr);CHKERRQ(ierr); 587 ierr = VecDestroy(&bl);CHKERRQ(ierr); 588 ierr = VecDestroy(&br);CHKERRQ(ierr); 589 PetscFunctionReturn(0); 590 } 591 592 #undef __FUNCT__ 593 #define __FUNCT__ "MatGetVecs_Nest" 594 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 595 { 596 Mat_Nest *bA = (Mat_Nest*)A->data; 597 Vec *L,*R; 598 MPI_Comm comm; 599 PetscInt i,j; 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 comm = ((PetscObject)A)->comm; 604 if (right) { 605 /* allocate R */ 606 ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr); 607 /* Create the right vectors */ 608 for (j=0; j<bA->nc; j++) { 609 for (i=0; i<bA->nr; i++) { 610 if (bA->m[i][j]) { 611 ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 612 break; 613 } 614 } 615 if (i==bA->nr) { 616 /* have an empty column */ 617 SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 618 } 619 } 620 ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 621 /* hand back control to the nest vector */ 622 for (j=0; j<bA->nc; j++) { 623 ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 624 } 625 ierr = PetscFree(R);CHKERRQ(ierr); 626 } 627 628 if (left) { 629 /* allocate L */ 630 ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr); 631 /* Create the left vectors */ 632 for (i=0; i<bA->nr; i++) { 633 for (j=0; j<bA->nc; j++) { 634 if (bA->m[i][j]) { 635 ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 636 break; 637 } 638 } 639 if (j==bA->nc) { 640 /* have an empty row */ 641 SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 642 } 643 } 644 645 ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 646 for (i=0; i<bA->nr; i++) { 647 ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 648 } 649 650 ierr = PetscFree(L);CHKERRQ(ierr); 651 } 652 PetscFunctionReturn(0); 653 } 654 655 #undef __FUNCT__ 656 #define __FUNCT__ "MatView_Nest" 657 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 658 { 659 Mat_Nest *bA = (Mat_Nest*)A->data; 660 PetscBool isascii; 661 PetscInt i,j; 662 PetscErrorCode ierr; 663 664 PetscFunctionBegin; 665 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 666 if (isascii) { 667 668 PetscViewerASCIIPrintf(viewer,"Matrix object: \n" ); 669 PetscViewerASCIIPushTab( viewer ); /* push0 */ 670 PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 671 672 PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" ); 673 for (i=0; i<bA->nr; i++) { 674 for (j=0; j<bA->nc; j++) { 675 const MatType type; 676 const char *name; 677 PetscInt NR,NC; 678 PetscBool isNest = PETSC_FALSE; 679 680 if (!bA->m[i][j]) { 681 PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j); 682 continue; 683 } 684 ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 685 ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr); 686 name = ((PetscObject)bA->m[i][j])->prefix; 687 ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 688 689 PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC); 690 691 if (isNest) { 692 PetscViewerASCIIPushTab( viewer ); /* push1 */ 693 ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr); 694 PetscViewerASCIIPopTab(viewer); /* pop1 */ 695 } 696 } 697 } 698 PetscViewerASCIIPopTab(viewer); /* pop0 */ 699 } 700 PetscFunctionReturn(0); 701 } 702 703 #undef __FUNCT__ 704 #define __FUNCT__ "MatZeroEntries_Nest" 705 static PetscErrorCode MatZeroEntries_Nest(Mat A) 706 { 707 Mat_Nest *bA = (Mat_Nest*)A->data; 708 PetscInt i,j; 709 PetscErrorCode ierr; 710 711 PetscFunctionBegin; 712 for (i=0; i<bA->nr; i++) { 713 for (j=0; j<bA->nc; j++) { 714 if (!bA->m[i][j]) continue; 715 ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 716 } 717 } 718 PetscFunctionReturn(0); 719 } 720 721 #undef __FUNCT__ 722 #define __FUNCT__ "MatDuplicate_Nest" 723 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 724 { 725 Mat_Nest *bA = (Mat_Nest*)A->data; 726 Mat *b; 727 PetscInt i,j,nr = bA->nr,nc = bA->nc; 728 PetscErrorCode ierr; 729 730 PetscFunctionBegin; 731 ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 732 for (i=0; i<nr; i++) { 733 for (j=0; j<nc; j++) { 734 if (bA->m[i][j]) { 735 ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 736 } else { 737 b[i*nc+j] = PETSC_NULL; 738 } 739 } 740 } 741 ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 742 /* Give the new MatNest exclusive ownership */ 743 for (i=0; i<nr*nc; i++) { 744 ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 745 } 746 ierr = PetscFree(b);CHKERRQ(ierr); 747 748 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 749 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 750 PetscFunctionReturn(0); 751 } 752 753 /* nest api */ 754 EXTERN_C_BEGIN 755 #undef __FUNCT__ 756 #define __FUNCT__ "MatNestGetSubMat_Nest" 757 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 758 { 759 Mat_Nest *bA = (Mat_Nest*)A->data; 760 PetscFunctionBegin; 761 if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 762 if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 763 *mat = bA->m[idxm][jdxm]; 764 PetscFunctionReturn(0); 765 } 766 EXTERN_C_END 767 768 #undef __FUNCT__ 769 #define __FUNCT__ "MatNestGetSubMat" 770 /*@C 771 MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 772 773 Not collective 774 775 Input Parameters: 776 + A - nest matrix 777 . idxm - index of the matrix within the nest matrix 778 - jdxm - index of the matrix within the nest matrix 779 780 Output Parameter: 781 . sub - matrix at index idxm,jdxm within the nest matrix 782 783 Level: developer 784 785 .seealso: MatNestGetSize(), MatNestGetSubMats() 786 @*/ 787 PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 788 { 789 PetscErrorCode ierr; 790 791 PetscFunctionBegin; 792 ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 793 PetscFunctionReturn(0); 794 } 795 796 EXTERN_C_BEGIN 797 #undef __FUNCT__ 798 #define __FUNCT__ "MatNestSetSubMat_Nest" 799 PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 800 { 801 Mat_Nest *bA = (Mat_Nest*)A->data; 802 PetscInt m,n,M,N,mi,ni,Mi,Ni; 803 PetscErrorCode ierr; 804 805 PetscFunctionBegin; 806 if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 807 if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 808 ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 809 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 810 ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 811 ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 812 ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 813 ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 814 if (M != Mi || N != Ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni); 815 if (m != mi || n != ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni); 816 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 817 ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 818 bA->m[idxm][jdxm] = mat; 819 PetscFunctionReturn(0); 820 } 821 EXTERN_C_END 822 823 #undef __FUNCT__ 824 #define __FUNCT__ "MatNestSetSubMat" 825 /*@C 826 MatNestSetSubMat - Set a single submatrix in the nest matrix. 827 828 Logically collective on the submatrix communicator 829 830 Input Parameters: 831 + A - nest matrix 832 . idxm - index of the matrix within the nest matrix 833 . jdxm - index of the matrix within the nest matrix 834 - sub - matrix at index idxm,jdxm within the nest matrix 835 836 Notes: 837 The new submatrix must have the same size and communicator as that block of the nest. 838 839 This increments the reference count of the submatrix. 840 841 Level: developer 842 843 .seealso: MatNestSetSubMats(), MatNestGetSubMat() 844 @*/ 845 PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 846 { 847 PetscErrorCode ierr; 848 849 PetscFunctionBegin; 850 ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 851 PetscFunctionReturn(0); 852 } 853 854 EXTERN_C_BEGIN 855 #undef __FUNCT__ 856 #define __FUNCT__ "MatNestGetSubMats_Nest" 857 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 858 { 859 Mat_Nest *bA = (Mat_Nest*)A->data; 860 PetscFunctionBegin; 861 if (M) { *M = bA->nr; } 862 if (N) { *N = bA->nc; } 863 if (mat) { *mat = bA->m; } 864 PetscFunctionReturn(0); 865 } 866 EXTERN_C_END 867 868 #undef __FUNCT__ 869 #define __FUNCT__ "MatNestGetSubMats" 870 /*@C 871 MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 872 873 Not collective 874 875 Input Parameters: 876 . A - nest matrix 877 878 Output Parameter: 879 + M - number of rows in the nest matrix 880 . N - number of cols in the nest matrix 881 - mat - 2d array of matrices 882 883 Notes: 884 885 The user should not free the array mat. 886 887 Level: developer 888 889 .seealso: MatNestGetSize(), MatNestGetSubMat() 890 @*/ 891 PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 892 { 893 PetscErrorCode ierr; 894 895 PetscFunctionBegin; 896 ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 897 PetscFunctionReturn(0); 898 } 899 900 EXTERN_C_BEGIN 901 #undef __FUNCT__ 902 #define __FUNCT__ "MatNestGetSize_Nest" 903 PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 904 { 905 Mat_Nest *bA = (Mat_Nest*)A->data; 906 907 PetscFunctionBegin; 908 if (M) { *M = bA->nr; } 909 if (N) { *N = bA->nc; } 910 PetscFunctionReturn(0); 911 } 912 EXTERN_C_END 913 914 #undef __FUNCT__ 915 #define __FUNCT__ "MatNestGetSize" 916 /*@C 917 MatNestGetSize - Returns the size of the nest matrix. 918 919 Not collective 920 921 Input Parameters: 922 . A - nest matrix 923 924 Output Parameter: 925 + M - number of rows in the nested mat 926 - N - number of cols in the nested mat 927 928 Notes: 929 930 Level: developer 931 932 .seealso: MatNestGetSubMat(), MatNestGetSubMats() 933 @*/ 934 PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 935 { 936 PetscErrorCode ierr; 937 938 PetscFunctionBegin; 939 ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 940 PetscFunctionReturn(0); 941 } 942 943 EXTERN_C_BEGIN 944 #undef __FUNCT__ 945 #define __FUNCT__ "MatNestSetVecType_Nest" 946 PetscErrorCode MatNestSetVecType_Nest(Mat A,const VecType vtype) 947 { 948 PetscErrorCode ierr; 949 PetscBool flg; 950 951 PetscFunctionBegin; 952 ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 953 /* In reality, this only distinguishes VECNEST and "other" */ 954 A->ops->getvecs = flg ? MatGetVecs_Nest : 0; 955 A->ops->getdiagonal = flg ? MatGetDiagonal_Nest : 0; 956 A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0; 957 958 PetscFunctionReturn(0); 959 } 960 EXTERN_C_END 961 962 #undef __FUNCT__ 963 #define __FUNCT__ "MatNestSetVecType" 964 /*@C 965 MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 966 967 Not collective 968 969 Input Parameters: 970 + A - nest matrix 971 - vtype - type to use for creating vectors 972 973 Notes: 974 975 Level: developer 976 977 .seealso: MatGetVecs() 978 @*/ 979 PetscErrorCode MatNestSetVecType(Mat A,const VecType vtype) 980 { 981 PetscErrorCode ierr; 982 983 PetscFunctionBegin; 984 ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr); 985 PetscFunctionReturn(0); 986 } 987 988 EXTERN_C_BEGIN 989 #undef __FUNCT__ 990 #define __FUNCT__ "MatNestSetSubMats_Nest" 991 PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 992 { 993 Mat_Nest *s = (Mat_Nest*)A->data; 994 PetscInt i,j,m,n,M,N; 995 PetscErrorCode ierr; 996 997 PetscFunctionBegin; 998 s->nr = nr; 999 s->nc = nc; 1000 1001 /* Create space for submatrices */ 1002 ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 1003 for (i=0; i<nr; i++) { 1004 ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 1005 } 1006 for (i=0; i<nr; i++) { 1007 for (j=0; j<nc; j++) { 1008 s->m[i][j] = a[i*nc+j]; 1009 if (a[i*nc+j]) { 1010 ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1011 } 1012 } 1013 } 1014 1015 ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1016 1017 ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 1018 ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 1019 for (i=0; i<nr; i++) s->row_len[i]=-1; 1020 for (j=0; j<nc; j++) s->col_len[j]=-1; 1021 1022 ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1023 1024 ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1025 ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1026 ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 1027 ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1028 ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1029 ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 1030 1031 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1032 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1033 1034 ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 1035 ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr); 1036 ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr); 1037 PetscFunctionReturn(0); 1038 } 1039 EXTERN_C_END 1040 1041 #undef __FUNCT__ 1042 #define __FUNCT__ "MatNestSetSubMats" 1043 /*@ 1044 MatNestSetSubMats - Sets the nested submatrices 1045 1046 Collective on Mat 1047 1048 Input Parameter: 1049 + N - nested matrix 1050 . nr - number of nested row blocks 1051 . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1052 . nc - number of nested column blocks 1053 . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1054 - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1055 1056 Level: advanced 1057 1058 .seealso: MatCreateNest(), MATNEST 1059 @*/ 1060 PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1061 { 1062 PetscErrorCode ierr; 1063 PetscInt i; 1064 1065 PetscFunctionBegin; 1066 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1067 if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1068 if (nr && is_row) { 1069 PetscValidPointer(is_row,3); 1070 for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1071 } 1072 if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 1073 if (nc && is_col) { 1074 PetscValidPointer(is_col,5); 1075 for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1076 } 1077 if (nr*nc) PetscValidPointer(a,6); 1078 ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr); 1079 PetscFunctionReturn(0); 1080 } 1081 1082 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1083 /* 1084 nprocessors = NP 1085 Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) ) 1086 proc 0: => (g_0,h_0,) 1087 proc 1: => (g_1,h_1,) 1088 ... 1089 proc nprocs-1: => (g_NP-1,h_NP-1,) 1090 1091 proc 0: proc 1: proc nprocs-1: 1092 is[0] = ( 0,1,2,...,nlocal(g_0)-1 ) ( 0,1,...,nlocal(g_1)-1 ) ( 0,1,...,nlocal(g_NP-1) ) 1093 1094 proc 0: 1095 is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 ) 1096 proc 1: 1097 is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 ) 1098 1099 proc NP-1: 1100 is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 ) 1101 */ 1102 #undef __FUNCT__ 1103 #define __FUNCT__ "MatSetUp_NestIS_Private" 1104 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1105 { 1106 Mat_Nest *vs = (Mat_Nest*)A->data; 1107 PetscInt i,j,offset,n,nsum,bs; 1108 PetscErrorCode ierr; 1109 Mat sub; 1110 1111 PetscFunctionBegin; 1112 ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 1113 ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1114 if (is_row) { /* valid IS is passed in */ 1115 /* refs on is[] are incremeneted */ 1116 for (i=0; i<vs->nr; i++) { 1117 ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1118 vs->isglobal.row[i] = is_row[i]; 1119 } 1120 } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 1121 nsum = 0; 1122 for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 1123 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1124 if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 1125 ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 1126 if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 1127 nsum += n; 1128 } 1129 offset = 0; 1130 ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1131 for (i=0; i<vs->nr; i++) { 1132 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1133 ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 1134 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1135 ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1136 ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 1137 offset += n; 1138 } 1139 } 1140 1141 if (is_col) { /* valid IS is passed in */ 1142 /* refs on is[] are incremeneted */ 1143 for (j=0; j<vs->nc; j++) { 1144 ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1145 vs->isglobal.col[j] = is_col[j]; 1146 } 1147 } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 1148 offset = A->cmap->rstart; 1149 nsum = 0; 1150 for (j=0; j<vs->nc; j++) { 1151 ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1152 if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 1153 ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 1154 if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 1155 nsum += n; 1156 } 1157 offset = 0; 1158 ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1159 for (j=0; j<vs->nc; j++) { 1160 ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1161 ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 1162 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1163 ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1164 ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 1165 offset += n; 1166 } 1167 } 1168 1169 /* Set up the local ISs */ 1170 ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1171 ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1172 for (i=0,offset=0; i<vs->nr; i++) { 1173 IS isloc; 1174 ISLocalToGlobalMapping rmap = PETSC_NULL; 1175 PetscInt nlocal,bs; 1176 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1177 if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);} 1178 if (rmap) { 1179 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1180 ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1181 ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1182 ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1183 } else { 1184 nlocal = 0; 1185 isloc = PETSC_NULL; 1186 } 1187 vs->islocal.row[i] = isloc; 1188 offset += nlocal; 1189 } 1190 for (i=0,offset=0; i<vs->nc; i++) { 1191 IS isloc; 1192 ISLocalToGlobalMapping cmap = PETSC_NULL; 1193 PetscInt nlocal,bs; 1194 ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 1195 if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);} 1196 if (cmap) { 1197 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1198 ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1199 ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1200 ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1201 } else { 1202 nlocal = 0; 1203 isloc = PETSC_NULL; 1204 } 1205 vs->islocal.col[i] = isloc; 1206 offset += nlocal; 1207 } 1208 1209 #if defined(PETSC_USE_DEBUG) 1210 for (i=0; i<vs->nr; i++) { 1211 for (j=0; j<vs->nc; j++) { 1212 PetscInt m,n,M,N,mi,ni,Mi,Ni; 1213 Mat B = vs->m[i][j]; 1214 if (!B) continue; 1215 ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 1216 ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 1217 ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 1218 ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 1219 ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 1220 ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 1221 if (M != Mi || N != Ni) SETERRQ6(((PetscObject)sub)->comm,PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni); 1222 if (m != mi || n != ni) SETERRQ6(((PetscObject)sub)->comm,PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni); 1223 } 1224 } 1225 #endif 1226 PetscFunctionReturn(0); 1227 } 1228 1229 #undef __FUNCT__ 1230 #define __FUNCT__ "MatCreateNest" 1231 /*@ 1232 MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1233 1234 Collective on Mat 1235 1236 Input Parameter: 1237 + comm - Communicator for the new Mat 1238 . nr - number of nested row blocks 1239 . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1240 . nc - number of nested column blocks 1241 . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1242 - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1243 1244 Output Parameter: 1245 . B - new matrix 1246 1247 Level: advanced 1248 1249 .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST 1250 @*/ 1251 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1252 { 1253 Mat A; 1254 PetscErrorCode ierr; 1255 1256 PetscFunctionBegin; 1257 *B = 0; 1258 ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1259 ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 1260 ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1261 *B = A; 1262 PetscFunctionReturn(0); 1263 } 1264 1265 /*MC 1266 MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1267 1268 Level: intermediate 1269 1270 Notes: 1271 This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1272 It allows the use of symmetric and block formats for parts of multi-physics simulations. 1273 It is usually used with DMComposite and DMGetMatrix() 1274 1275 .seealso: MatCreate(), MatType, MatCreateNest() 1276 M*/ 1277 EXTERN_C_BEGIN 1278 #undef __FUNCT__ 1279 #define __FUNCT__ "MatCreate_Nest" 1280 PetscErrorCode MatCreate_Nest(Mat A) 1281 { 1282 Mat_Nest *s; 1283 PetscErrorCode ierr; 1284 1285 PetscFunctionBegin; 1286 ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr); 1287 A->data = (void*)s; 1288 s->nr = s->nc = -1; 1289 s->m = PETSC_NULL; 1290 1291 ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 1292 A->ops->mult = MatMult_Nest; 1293 A->ops->multadd = 0; 1294 A->ops->multtranspose = MatMultTranspose_Nest; 1295 A->ops->multtransposeadd = 0; 1296 A->ops->assemblybegin = MatAssemblyBegin_Nest; 1297 A->ops->assemblyend = MatAssemblyEnd_Nest; 1298 A->ops->zeroentries = MatZeroEntries_Nest; 1299 A->ops->duplicate = MatDuplicate_Nest; 1300 A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1301 A->ops->destroy = MatDestroy_Nest; 1302 A->ops->view = MatView_Nest; 1303 A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1304 A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1305 A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1306 A->ops->getdiagonal = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */ 1307 A->ops->diagonalscale = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */ 1308 1309 A->spptr = 0; 1310 A->same_nonzero = PETSC_FALSE; 1311 A->assembled = PETSC_FALSE; 1312 1313 /* expose Nest api's */ 1314 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "MatNestGetSubMat_Nest", MatNestGetSubMat_Nest);CHKERRQ(ierr); 1315 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "MatNestSetSubMat_Nest", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1316 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "MatNestGetSubMats_Nest", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1317 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "MatNestGetSize_Nest", MatNestGetSize_Nest);CHKERRQ(ierr); 1318 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "MatNestSetVecType_Nest", MatNestSetVecType_Nest);CHKERRQ(ierr); 1319 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "MatNestSetSubMats_Nest", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1320 1321 ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1322 PetscFunctionReturn(0); 1323 } 1324 EXTERN_C_END 1325