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