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