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