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