1 #define PETSCMAT_DLL 2 3 #include "matnestimpl.h" /*I "petscmat.h" I*/ 4 5 /* private functions */ 6 #undef __FUNCT__ 7 #define __FUNCT__ "MatNestGetSize_Recursive" 8 static PetscErrorCode MatNestGetSize_Recursive(Mat A,PetscBool globalsize,PetscInt *M,PetscInt *N) 9 { 10 Mat_Nest *bA = (Mat_Nest*)A->data; 11 PetscInt sizeM,sizeN,i,j,index = -1; 12 PetscErrorCode ierr; 13 14 PetscFunctionBegin; 15 /* rows */ 16 /* now descend recursively */ 17 for (i=0; i<bA->nr; i++) { 18 for (j=0; j<bA->nc; j++) { 19 if (bA->m[i][j]) { index = j; break; } 20 } 21 if (bA->m[i][index]) { 22 if (globalsize) { ierr = MatGetSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); } 23 else { ierr = MatGetLocalSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); } 24 *M = *M + sizeM; 25 } 26 } 27 28 /* cols */ 29 /* now descend recursively */ 30 for (j=0; j<bA->nc; j++) { 31 for (i=0; i<bA->nr; i++) { 32 if (bA->m[i][j]) { index = i; break; } 33 } 34 if (bA->m[index][j]) { 35 if (globalsize) { ierr = MatGetSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); } 36 else { ierr = MatGetLocalSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); } 37 *N = *N + sizeN; 38 } 39 } 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2" 45 PETSC_UNUSED 46 static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y) 47 { 48 PetscBool isANest, isxNest, isyNest; 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 isANest = isxNest = isyNest = PETSC_FALSE; 53 ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr); 54 ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr); 55 ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr); 56 57 if (isANest && isxNest && isyNest) { 58 PetscFunctionReturn(0); 59 } else { 60 SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations"); 61 PetscFunctionReturn(0); 62 } 63 PetscFunctionReturn(0); 64 } 65 66 /* 67 This function is called every time we insert a sub matrix the Nest. 68 It ensures that every Nest along row r, and col c has the same dimensions 69 as the submat being inserted. 70 */ 71 #undef __FUNCT__ 72 #define __FUNCT__ "PETSc_MatNest_CheckConsistency" 73 PETSC_UNUSED 74 static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c) 75 { 76 Mat_Nest *b = (Mat_Nest*)A->data; 77 PetscInt i,j; 78 PetscInt nr,nc; 79 PetscInt sM,sN,M,N; 80 Mat mat; 81 PetscErrorCode ierr; 82 83 PetscFunctionBegin; 84 nr = b->nr; 85 nc = b->nc; 86 ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr); 87 for (i=0; i<nr; i++) { 88 mat = b->m[i][c]; 89 if (mat) { 90 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 91 /* Check columns match */ 92 if (sN != N) { 93 SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N ); 94 } 95 } 96 } 97 98 for (j=0; j<nc; j++) { 99 mat = b->m[r][j]; 100 if (mat) { 101 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 102 /* Check rows match */ 103 if (sM != M) { 104 SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M ); 105 } 106 } 107 } 108 PetscFunctionReturn(0); 109 } 110 111 /* 112 Checks if the row/col's contain a complete set of IS's. 113 Once they are filled, the offsets are computed. This saves having to update 114 the index set entries each time we insert something new. 115 */ 116 #undef __FUNCT__ 117 #define __FUNCT__ "PETSc_MatNest_UpdateStructure" 118 PETSC_UNUSED 119 static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx) 120 { 121 Mat_Nest *b = (Mat_Nest*)A->data; 122 PetscInt m,n,i,j; 123 PetscBool fullrow,fullcol; 124 PetscErrorCode ierr; 125 126 PetscFunctionBegin; 127 ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr); 128 b->row_len[ridx] = m; 129 b->col_len[cidx] = n; 130 131 fullrow = PETSC_TRUE; 132 for (i=0; i<b->nr; i++) { 133 if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; } 134 } 135 if ( (fullrow) && (!b->isglobal.row[0]) ){ 136 PetscInt cnt = 0; 137 for (i=0; i<b->nr; i++) { 138 ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr); 139 cnt = cnt + b->row_len[i]; 140 } 141 } 142 143 fullcol = PETSC_TRUE; 144 for (j=0; j<b->nc; j++) { 145 if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; } 146 } 147 if( (fullcol) && (!b->isglobal.col[0]) ){ 148 PetscInt cnt = 0; 149 for (j=0; j<b->nc; j++) { 150 ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr); 151 cnt = cnt + b->col_len[j]; 152 } 153 } 154 PetscFunctionReturn(0); 155 } 156 157 /* operations */ 158 #undef __FUNCT__ 159 #define __FUNCT__ "MatMult_Nest" 160 static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 161 { 162 Mat_Nest *bA = (Mat_Nest*)A->data; 163 Vec *bx = bA->right,*by = bA->left; 164 PetscInt i,j,nr = bA->nr,nc = bA->nc; 165 PetscErrorCode ierr; 166 167 PetscFunctionBegin; 168 for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 169 for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 170 for (i=0; i<nr; i++) { 171 ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 172 for (j=0; j<nc; j++) { 173 if (!bA->m[i][j]) continue; 174 /* y[i] <- y[i] + A[i][j] * x[j] */ 175 ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 176 } 177 } 178 for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 179 for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 180 PetscFunctionReturn(0); 181 } 182 183 #undef __FUNCT__ 184 #define __FUNCT__ "MatMultTranspose_Nest" 185 static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 186 { 187 Mat_Nest *bA = (Mat_Nest*)A->data; 188 Vec *bx = bA->left,*by = bA->right; 189 PetscInt i,j,nr = bA->nr,nc = bA->nc; 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 if (A->symmetric) { 194 ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr); 195 PetscFunctionReturn(0); 196 } 197 for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 198 for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 199 for (i=0; i<nr; i++) { 200 ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 201 for (j=0; j<nc; j++) { 202 if (!bA->m[j][i]) continue; 203 /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */ 204 ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr); 205 } 206 } 207 for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 208 for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 209 PetscFunctionReturn(0); 210 } 211 212 #undef __FUNCT__ 213 #define __FUNCT__ "MatNestDestroyISList" 214 static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 215 { 216 PetscErrorCode ierr; 217 IS *lst = *list; 218 PetscInt i; 219 220 PetscFunctionBegin; 221 if (!lst) PetscFunctionReturn(0); 222 for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(lst[i]);CHKERRQ(ierr);} 223 ierr = PetscFree(lst);CHKERRQ(ierr); 224 *list = PETSC_NULL; 225 PetscFunctionReturn(0); 226 } 227 228 #undef __FUNCT__ 229 #define __FUNCT__ "MatDestroy_Nest" 230 static PetscErrorCode MatDestroy_Nest(Mat A) 231 { 232 Mat_Nest *vs = (Mat_Nest*)A->data; 233 PetscInt i,j; 234 PetscErrorCode ierr; 235 236 PetscFunctionBegin; 237 /* release the matrices and the place holders */ 238 ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 239 ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 240 ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 241 ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 242 243 ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 244 ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 245 246 ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 247 248 /* release the matrices and the place holders */ 249 if (vs->m) { 250 for (i=0; i<vs->nr; i++) { 251 for (j=0; j<vs->nc; j++) { 252 253 if (vs->m[i][j]) { 254 ierr = MatDestroy(vs->m[i][j]);CHKERRQ(ierr); 255 vs->m[i][j] = PETSC_NULL; 256 } 257 } 258 ierr = PetscFree( vs->m[i] );CHKERRQ(ierr); 259 vs->m[i] = PETSC_NULL; 260 } 261 ierr = PetscFree(vs->m);CHKERRQ(ierr); 262 vs->m = PETSC_NULL; 263 } 264 ierr = PetscFree(vs);CHKERRQ(ierr); 265 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 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 *B = PETSC_NULL; 435 PetscFunctionReturn(0); 436 } 437 438 #undef __FUNCT__ 439 #define __FUNCT__ "MatGetVecs_Nest" 440 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 441 { 442 Mat_Nest *bA = (Mat_Nest*)A->data; 443 Vec *L,*R; 444 MPI_Comm comm; 445 PetscInt i,j; 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 comm = ((PetscObject)A)->comm; 450 if (right) { 451 /* allocate R */ 452 ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr); 453 /* Create the right vectors */ 454 for (j=0; j<bA->nc; j++) { 455 for (i=0; i<bA->nr; i++) { 456 if (bA->m[i][j]) { 457 ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 458 break; 459 } 460 } 461 if (i==bA->nr) { 462 /* have an empty column */ 463 SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 464 } 465 } 466 ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 467 /* hand back control to the nest vector */ 468 for (j=0; j<bA->nc; j++) { 469 ierr = VecDestroy(R[j]);CHKERRQ(ierr); 470 } 471 ierr = PetscFree(R);CHKERRQ(ierr); 472 } 473 474 if (left) { 475 /* allocate L */ 476 ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr); 477 /* Create the left vectors */ 478 for (i=0; i<bA->nr; i++) { 479 for (j=0; j<bA->nc; j++) { 480 if (bA->m[i][j]) { 481 ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 482 break; 483 } 484 } 485 if (j==bA->nc) { 486 /* have an empty row */ 487 SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 488 } 489 } 490 491 ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 492 for (i=0; i<bA->nr; i++) { 493 ierr = VecDestroy(L[i]);CHKERRQ(ierr); 494 } 495 496 ierr = PetscFree(L);CHKERRQ(ierr); 497 } 498 PetscFunctionReturn(0); 499 } 500 501 #undef __FUNCT__ 502 #define __FUNCT__ "MatView_Nest" 503 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 504 { 505 Mat_Nest *bA = (Mat_Nest*)A->data; 506 PetscBool isascii; 507 PetscInt i,j; 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 512 if (isascii) { 513 514 PetscViewerASCIIPrintf(viewer,"Matrix object: \n" ); 515 PetscViewerASCIIPushTab( viewer ); /* push0 */ 516 PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 517 518 PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" ); 519 for (i=0; i<bA->nr; i++) { 520 for (j=0; j<bA->nc; j++) { 521 const MatType type; 522 const char *name; 523 PetscInt NR,NC; 524 PetscBool isNest = PETSC_FALSE; 525 526 if (!bA->m[i][j]) { 527 PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j); 528 continue; 529 } 530 ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 531 ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr); 532 name = ((PetscObject)bA->m[i][j])->prefix; 533 ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 534 535 PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC); 536 537 if (isNest) { 538 PetscViewerASCIIPushTab( viewer ); /* push1 */ 539 ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr); 540 PetscViewerASCIIPopTab(viewer); /* pop1 */ 541 } 542 } 543 } 544 PetscViewerASCIIPopTab(viewer); /* pop0 */ 545 } 546 PetscFunctionReturn(0); 547 } 548 549 #undef __FUNCT__ 550 #define __FUNCT__ "MatZeroEntries_Nest" 551 static PetscErrorCode MatZeroEntries_Nest(Mat A) 552 { 553 Mat_Nest *bA = (Mat_Nest*)A->data; 554 PetscInt i,j; 555 PetscErrorCode ierr; 556 557 PetscFunctionBegin; 558 for (i=0; i<bA->nr; i++) { 559 for (j=0; j<bA->nc; j++) { 560 if (!bA->m[i][j]) continue; 561 ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 562 } 563 } 564 PetscFunctionReturn(0); 565 } 566 567 #undef __FUNCT__ 568 #define __FUNCT__ "MatDuplicate_Nest" 569 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 570 { 571 Mat_Nest *bA = (Mat_Nest*)A->data; 572 Mat *b; 573 PetscInt i,j,nr = bA->nr,nc = bA->nc; 574 PetscErrorCode ierr; 575 576 PetscFunctionBegin; 577 ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 578 for (i=0; i<nr; i++) { 579 for (j=0; j<nc; j++) { 580 if (bA->m[i][j]) { 581 ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 582 } else { 583 b[i*nc+j] = PETSC_NULL; 584 } 585 } 586 } 587 ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 588 /* Give the new MatNest exclusive ownership */ 589 for (i=0; i<nr*nc; i++) { 590 if (b[i]) {ierr = MatDestroy(b[i]);CHKERRQ(ierr);} 591 } 592 ierr = PetscFree(b);CHKERRQ(ierr); 593 594 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 595 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 596 PetscFunctionReturn(0); 597 } 598 599 /* nest api */ 600 EXTERN_C_BEGIN 601 #undef __FUNCT__ 602 #define __FUNCT__ "MatNestGetSubMat_Nest" 603 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 604 { 605 Mat_Nest *bA = (Mat_Nest*)A->data; 606 PetscFunctionBegin; 607 if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 608 if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 609 *mat = bA->m[idxm][jdxm]; 610 PetscFunctionReturn(0); 611 } 612 EXTERN_C_END 613 614 #undef __FUNCT__ 615 #define __FUNCT__ "MatNestGetSubMat" 616 /*@C 617 MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 618 619 Not collective 620 621 Input Parameters: 622 . A - nest matrix 623 . idxm - index of the matrix within the nest matrix 624 . jdxm - index of the matrix within the nest matrix 625 626 Output Parameter: 627 . sub - matrix at index idxm,jdxm within the nest matrix 628 629 Notes: 630 631 Level: developer 632 633 .seealso: MatNestGetSize(), MatNestGetSubMats() 634 @*/ 635 PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 636 { 637 PetscErrorCode ierr; 638 639 PetscFunctionBegin; 640 ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 641 PetscFunctionReturn(0); 642 } 643 644 EXTERN_C_BEGIN 645 #undef __FUNCT__ 646 #define __FUNCT__ "MatNestGetSubMats_Nest" 647 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 648 { 649 Mat_Nest *bA = (Mat_Nest*)A->data; 650 PetscFunctionBegin; 651 if (M) { *M = bA->nr; } 652 if (N) { *N = bA->nc; } 653 if (mat) { *mat = bA->m; } 654 PetscFunctionReturn(0); 655 } 656 EXTERN_C_END 657 658 #undef __FUNCT__ 659 #define __FUNCT__ "MatNestGetSubMats" 660 /*@C 661 MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 662 663 Not collective 664 665 Input Parameters: 666 . A - nest matri 667 668 Output Parameter: 669 . M - number of rows in the nest matrix 670 . N - number of cols in the nest matrix 671 . mat - 2d array of matrices 672 673 Notes: 674 675 The user should not free the array mat. 676 677 Level: developer 678 679 .seealso: MatNestGetSize(), MatNestGetSubMat() 680 @*/ 681 PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 682 { 683 PetscErrorCode ierr; 684 685 PetscFunctionBegin; 686 ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 687 PetscFunctionReturn(0); 688 } 689 690 EXTERN_C_BEGIN 691 #undef __FUNCT__ 692 #define __FUNCT__ "MatNestGetSize_Nest" 693 PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 694 { 695 Mat_Nest *bA = (Mat_Nest*)A->data; 696 697 PetscFunctionBegin; 698 if (M) { *M = bA->nr; } 699 if (N) { *N = bA->nc; } 700 PetscFunctionReturn(0); 701 } 702 EXTERN_C_END 703 704 #undef __FUNCT__ 705 #define __FUNCT__ "MatNestGetSize" 706 /*@C 707 MatNestGetSize - Returns the size of the nest matrix. 708 709 Not collective 710 711 Input Parameters: 712 . A - nest matrix 713 714 Output Parameter: 715 . M - number of rows in the nested mat 716 . N - number of cols in the nested mat 717 718 Notes: 719 720 Level: developer 721 722 .seealso: MatNestGetSubMat(), MatNestGetSubMats() 723 @*/ 724 PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 725 { 726 PetscErrorCode ierr; 727 728 PetscFunctionBegin; 729 ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 730 PetscFunctionReturn(0); 731 } 732 733 EXTERN_C_BEGIN 734 #undef __FUNCT__ 735 #define __FUNCT__ "MatNestSetVecType_Nest" 736 PetscErrorCode MatNestSetVecType_Nest(Mat A,const VecType vtype) 737 { 738 PetscErrorCode ierr; 739 PetscBool flg; 740 741 PetscFunctionBegin; 742 ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 743 /* In reality, this only distinguishes VECNEST and "other" */ 744 A->ops->getvecs = flg ? MatGetVecs_Nest : 0; 745 PetscFunctionReturn(0); 746 } 747 EXTERN_C_END 748 749 #undef __FUNCT__ 750 #define __FUNCT__ "MatNestSetVecType" 751 /*@C 752 MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 753 754 Not collective 755 756 Input Parameters: 757 + A - nest matrix 758 - vtype - type to use for creating vectors 759 760 Notes: 761 762 Level: developer 763 764 .seealso: MatGetVecs() 765 @*/ 766 PetscErrorCode MatNestSetVecType(Mat A,const VecType vtype) 767 { 768 PetscErrorCode ierr; 769 770 PetscFunctionBegin; 771 ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr); 772 PetscFunctionReturn(0); 773 } 774 775 /* constructor */ 776 #undef __FUNCT__ 777 #define __FUNCT__ "MatNestSetOps_Private" 778 static PetscErrorCode MatNestSetOps_Private(struct _MatOps* ops) 779 { 780 PetscFunctionBegin; 781 /* 0*/ 782 ops->setvalues = 0; 783 ops->getrow = 0; 784 ops->restorerow = 0; 785 ops->mult = MatMult_Nest; 786 ops->multadd = 0; 787 /* 5*/ 788 ops->multtranspose = MatMultTranspose_Nest; 789 ops->multtransposeadd = 0; 790 ops->solve = 0; 791 ops->solveadd = 0; 792 ops->solvetranspose = 0; 793 /*10*/ 794 ops->solvetransposeadd = 0; 795 ops->lufactor = 0; 796 ops->choleskyfactor = 0; 797 ops->sor = 0; 798 ops->transpose = 0; 799 /*15*/ 800 ops->getinfo = 0; 801 ops->equal = 0; 802 ops->getdiagonal = 0; 803 ops->diagonalscale = 0; 804 ops->norm = 0; 805 /*20*/ 806 ops->assemblybegin = MatAssemblyBegin_Nest; 807 ops->assemblyend = MatAssemblyEnd_Nest; 808 ops->setoption = 0; 809 ops->zeroentries = MatZeroEntries_Nest; 810 /*24*/ 811 ops->zerorows = 0; 812 ops->lufactorsymbolic = 0; 813 ops->lufactornumeric = 0; 814 ops->choleskyfactorsymbolic = 0; 815 ops->choleskyfactornumeric = 0; 816 /*29*/ 817 ops->setuppreallocation = 0; 818 ops->ilufactorsymbolic = 0; 819 ops->iccfactorsymbolic = 0; 820 ops->getarray = 0; 821 ops->restorearray = 0; 822 /*34*/ 823 ops->duplicate = MatDuplicate_Nest; 824 ops->forwardsolve = 0; 825 ops->backwardsolve = 0; 826 ops->ilufactor = 0; 827 ops->iccfactor = 0; 828 /*39*/ 829 ops->axpy = 0; 830 ops->getsubmatrices = 0; 831 ops->increaseoverlap = 0; 832 ops->getvalues = 0; 833 ops->copy = 0; 834 /*44*/ 835 ops->getrowmax = 0; 836 ops->scale = 0; 837 ops->shift = 0; 838 ops->diagonalset = 0; 839 ops->zerorowscolumns = 0; 840 /*49*/ 841 ops->setblocksize = 0; 842 ops->getrowij = 0; 843 ops->restorerowij = 0; 844 ops->getcolumnij = 0; 845 ops->restorecolumnij = 0; 846 /*54*/ 847 ops->fdcoloringcreate = 0; 848 ops->coloringpatch = 0; 849 ops->setunfactored = 0; 850 ops->permute = 0; 851 ops->setvaluesblocked = 0; 852 /*59*/ 853 ops->getsubmatrix = MatGetSubMatrix_Nest; 854 ops->destroy = MatDestroy_Nest; 855 ops->view = MatView_Nest; 856 ops->convertfrom = 0; 857 ops->usescaledform = 0; 858 /*64*/ 859 ops->scalesystem = 0; 860 ops->unscalesystem = 0; 861 ops->setlocaltoglobalmapping = 0; 862 ops->setvalueslocal = 0; 863 ops->zerorowslocal = 0; 864 /*69*/ 865 ops->getrowmaxabs = 0; 866 ops->getrowminabs = 0; 867 ops->convert = 0; 868 ops->setcoloring = 0; 869 ops->setvaluesadic = 0; 870 /* 74 */ 871 ops->setvaluesadifor = 0; 872 ops->fdcoloringapply = 0; 873 ops->setfromoptions = 0; 874 ops->multconstrained = 0; 875 ops->multtransposeconstrained = 0; 876 /*79*/ 877 ops->dummy = 0; 878 ops->mults = 0; 879 ops->solves = 0; 880 ops->getinertia = 0; 881 ops->load = 0; 882 /*84*/ 883 ops->issymmetric = 0; 884 ops->ishermitian = 0; 885 ops->isstructurallysymmetric = 0; 886 ops->dummy4 = 0; 887 ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 888 /*89*/ 889 ops->matmult = 0; 890 ops->matmultsymbolic = 0; 891 ops->matmultnumeric = 0; 892 ops->ptap = 0; 893 ops->ptapsymbolic = 0; 894 /*94*/ 895 ops->ptapnumeric = 0; 896 ops->matmulttranspose = 0; 897 ops->matmulttransposesymbolic = 0; 898 ops->matmulttransposenumeric = 0; 899 ops->ptapsymbolic_seqaij = 0; 900 /*99*/ 901 ops->ptapnumeric_seqaij = 0; 902 ops->ptapsymbolic_mpiaij = 0; 903 ops->ptapnumeric_mpiaij = 0; 904 ops->conjugate = 0; 905 ops->setsizes = 0; 906 /*104*/ 907 ops->setvaluesrow = 0; 908 ops->realpart = 0; 909 ops->imaginarypart = 0; 910 ops->getrowuppertriangular = 0; 911 ops->restorerowuppertriangular = 0; 912 /*109*/ 913 ops->matsolve = 0; 914 ops->getredundantmatrix = 0; 915 ops->getrowmin = 0; 916 ops->getcolumnvector = 0; 917 ops->missingdiagonal = 0; 918 /* 114 */ 919 ops->getseqnonzerostructure = 0; 920 ops->create = 0; 921 ops->getghosts = 0; 922 ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 923 ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 924 /* 119 */ 925 ops->multdiagonalblock = 0; 926 ops->hermitiantranspose = 0; 927 ops->multhermitiantranspose = 0; 928 ops->multhermitiantransposeadd = 0; 929 ops->getmultiprocblock = 0; 930 /* 124 */ 931 ops->dummy1 = 0; 932 ops->dummy2 = 0; 933 ops->dummy3 = 0; 934 ops->dummy4 = 0; 935 ops->getsubmatricesparallel = 0; 936 937 PetscFunctionReturn(0); 938 } 939 940 #undef __FUNCT__ 941 #define __FUNCT__ "MatSetUp_Nest_Private" 942 static PetscErrorCode MatSetUp_Nest_Private(Mat A,PetscInt nr,PetscInt nc,const Mat *sub) 943 { 944 Mat_Nest *ctx = (Mat_Nest*)A->data; 945 PetscInt i,j; 946 PetscErrorCode ierr; 947 948 PetscFunctionBegin; 949 if(ctx->setup_called) PetscFunctionReturn(0); 950 951 ctx->nr = nr; 952 ctx->nc = nc; 953 954 /* Create space */ 955 ierr = PetscMalloc(sizeof(Mat*)*ctx->nr,&ctx->m);CHKERRQ(ierr); 956 for (i=0; i<ctx->nr; i++) { 957 ierr = PetscMalloc(sizeof(Mat)*ctx->nc,&ctx->m[i]);CHKERRQ(ierr); 958 } 959 for (i=0; i<ctx->nr; i++) { 960 for (j=0; j<ctx->nc; j++) { 961 ctx->m[i][j] = sub[i*nc+j]; 962 if (sub[i*nc+j]) { 963 ierr = PetscObjectReference((PetscObject)sub[i*nc+j]);CHKERRQ(ierr); 964 } 965 } 966 } 967 968 ierr = PetscMalloc(sizeof(IS)*ctx->nr,&ctx->isglobal.row);CHKERRQ(ierr); 969 ierr = PetscMalloc(sizeof(IS)*ctx->nc,&ctx->isglobal.col);CHKERRQ(ierr); 970 971 ierr = PetscMalloc(sizeof(PetscInt)*ctx->nr,&ctx->row_len);CHKERRQ(ierr); 972 ierr = PetscMalloc(sizeof(PetscInt)*ctx->nc,&ctx->col_len);CHKERRQ(ierr); 973 for (i=0;i<ctx->nr;i++) ctx->row_len[i]=-1; 974 for (j=0;j<ctx->nc;j++) ctx->col_len[j]=-1; 975 976 ctx->setup_called = PETSC_TRUE; 977 978 PetscFunctionReturn(0); 979 } 980 981 982 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 983 /* 984 nprocessors = NP 985 Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) ) 986 proc 0: => (g_0,h_0,) 987 proc 1: => (g_1,h_1,) 988 ... 989 proc nprocs-1: => (g_NP-1,h_NP-1,) 990 991 proc 0: proc 1: proc nprocs-1: 992 is[0] = ( 0,1,2,...,nlocal(g_0)-1 ) ( 0,1,...,nlocal(g_1)-1 ) ( 0,1,...,nlocal(g_NP-1) ) 993 994 proc 0: 995 is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 ) 996 proc 1: 997 is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 ) 998 999 proc NP-1: 1000 is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 ) 1001 */ 1002 #undef __FUNCT__ 1003 #define __FUNCT__ "MatSetUp_NestIS_Private" 1004 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1005 { 1006 Mat_Nest *vs = (Mat_Nest*)A->data; 1007 PetscInt i,j,offset,n,bs; 1008 PetscErrorCode ierr; 1009 Mat sub; 1010 1011 PetscFunctionBegin; 1012 if (is_row) { /* valid IS is passed in */ 1013 /* refs on is[] are incremeneted */ 1014 for (i=0; i<vs->nr; i++) { 1015 ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1016 vs->isglobal.row[i] = is_row[i]; 1017 } 1018 } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 1019 offset = A->rmap->rstart; 1020 for (i=0; i<vs->nr; i++) { 1021 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1022 ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 1023 if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 1024 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1025 ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1026 ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 1027 offset += n; 1028 } 1029 } 1030 1031 if (is_col) { /* valid IS is passed in */ 1032 /* refs on is[] are incremeneted */ 1033 for (j=0; j<vs->nc; j++) { 1034 ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1035 vs->isglobal.col[j] = is_col[j]; 1036 } 1037 } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 1038 offset = A->cmap->rstart; 1039 for (j=0; j<vs->nc; j++) { 1040 ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1041 ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 1042 if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 1043 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1044 ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1045 ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 1046 offset += n; 1047 } 1048 } 1049 1050 /* Set up the local ISs */ 1051 ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1052 ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1053 for (i=0,offset=0; i<vs->nr; i++) { 1054 IS isloc; 1055 ISLocalToGlobalMapping rmap; 1056 PetscInt nlocal,bs; 1057 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1058 ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr); 1059 if (rmap) { 1060 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1061 ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1062 ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1063 ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1064 } else { 1065 nlocal = 0; 1066 isloc = PETSC_NULL; 1067 } 1068 vs->islocal.row[i] = isloc; 1069 offset += nlocal; 1070 } 1071 for (i=0,offset=0; i<vs->nr; i++) { 1072 IS isloc; 1073 ISLocalToGlobalMapping cmap; 1074 PetscInt nlocal,bs; 1075 ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 1076 ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr); 1077 if (cmap) { 1078 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1079 ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1080 ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1081 ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1082 } else { 1083 nlocal = 0; 1084 isloc = PETSC_NULL; 1085 } 1086 vs->islocal.col[i] = isloc; 1087 offset += nlocal; 1088 } 1089 PetscFunctionReturn(0); 1090 } 1091 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "MatCreateNest" 1094 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1095 { 1096 Mat A; 1097 Mat_Nest *s; 1098 PetscInt i,m,n,M,N; 1099 PetscErrorCode ierr; 1100 1101 PetscFunctionBegin; 1102 if (nr < 0) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1103 if (nr && is_row) { 1104 PetscValidPointer(is_row,3); 1105 for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1106 } 1107 if (nc < 0) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 1108 if (nc && is_row) { 1109 PetscValidPointer(is_col,5); 1110 for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1111 } 1112 if (nr*nc) PetscValidPointer(a,6); 1113 PetscValidPointer(B,7); 1114 1115 ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1116 1117 ierr = PetscMalloc( sizeof(Mat_Nest), &s );CHKERRQ(ierr); 1118 ierr = PetscMemzero( s, sizeof(Mat_Nest) );CHKERRQ(ierr); 1119 A->data = (void*)s; 1120 s->setup_called = PETSC_FALSE; 1121 s->nr = s->nc = -1; 1122 s->m = PETSC_NULL; 1123 1124 /* define the operations */ 1125 ierr = MatNestSetOps_Private(A->ops);CHKERRQ(ierr); 1126 A->spptr = 0; 1127 A->same_nonzero = PETSC_FALSE; 1128 A->assembled = PETSC_FALSE; 1129 1130 ierr = PetscObjectChangeTypeName((PetscObject) A, MATNEST );CHKERRQ(ierr); 1131 ierr = MatSetUp_Nest_Private(A,nr,nc,a);CHKERRQ(ierr); 1132 1133 m = n = 0; 1134 M = N = 0; 1135 ierr = MatNestGetSize_Recursive(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr); 1136 ierr = MatNestGetSize_Recursive(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr); 1137 1138 ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1139 ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1140 ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 1141 ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1142 ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1143 ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 1144 1145 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1146 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1147 1148 ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1149 ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 1150 1151 /* expose Nest api's */ 1152 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C","MatNestGetSubMat_Nest",MatNestGetSubMat_Nest);CHKERRQ(ierr); 1153 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C","MatNestGetSubMats_Nest",MatNestGetSubMats_Nest);CHKERRQ(ierr); 1154 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C","MatNestGetSize_Nest",MatNestGetSize_Nest);CHKERRQ(ierr); 1155 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C","MatNestSetVecType_Nest",MatNestSetVecType_Nest);CHKERRQ(ierr); 1156 1157 *B = A; 1158 PetscFunctionReturn(0); 1159 } 1160