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