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 = PetscMalloc1(bA->nc,&br);CHKERRQ(ierr); 469 if (r) { 470 for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 471 } else { 472 PetscMemzero(br, bA->nc*sizeof(Vec)); 473 } 474 bl = NULL; 475 for (i=0; i<bA->nr; i++) { 476 if (l) { 477 ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 478 } 479 for (j=0; j<bA->nc; j++) { 480 if (bA->m[i][j]) { 481 ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 482 } 483 } 484 if (l) { 485 ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 486 } 487 } 488 if (r) { 489 for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 490 } 491 ierr = PetscFree(br);CHKERRQ(ierr); 492 PetscFunctionReturn(0); 493 } 494 495 #undef __FUNCT__ 496 #define __FUNCT__ "MatScale_Nest" 497 static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 498 { 499 Mat_Nest *bA = (Mat_Nest*)A->data; 500 PetscInt i,j; 501 PetscErrorCode ierr; 502 503 PetscFunctionBegin; 504 for (i=0; i<bA->nr; i++) { 505 for (j=0; j<bA->nc; j++) { 506 if (bA->m[i][j]) { 507 ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 508 } 509 } 510 } 511 PetscFunctionReturn(0); 512 } 513 514 #undef __FUNCT__ 515 #define __FUNCT__ "MatShift_Nest" 516 static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 517 { 518 Mat_Nest *bA = (Mat_Nest*)A->data; 519 PetscInt i; 520 PetscErrorCode ierr; 521 522 PetscFunctionBegin; 523 for (i=0; i<bA->nr; i++) { 524 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); 525 ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 526 } 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "MatGetVecs_Nest" 532 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 533 { 534 Mat_Nest *bA = (Mat_Nest*)A->data; 535 Vec *L,*R; 536 MPI_Comm comm; 537 PetscInt i,j; 538 PetscErrorCode ierr; 539 540 PetscFunctionBegin; 541 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 542 if (right) { 543 /* allocate R */ 544 ierr = PetscMalloc(sizeof(Vec) * bA->nc, &R);CHKERRQ(ierr); 545 /* Create the right vectors */ 546 for (j=0; j<bA->nc; j++) { 547 for (i=0; i<bA->nr; i++) { 548 if (bA->m[i][j]) { 549 ierr = MatGetVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr); 550 break; 551 } 552 } 553 if (i==bA->nr) { 554 /* have an empty column */ 555 SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 556 } 557 } 558 ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 559 /* hand back control to the nest vector */ 560 for (j=0; j<bA->nc; j++) { 561 ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 562 } 563 ierr = PetscFree(R);CHKERRQ(ierr); 564 } 565 566 if (left) { 567 /* allocate L */ 568 ierr = PetscMalloc(sizeof(Vec) * bA->nr, &L);CHKERRQ(ierr); 569 /* Create the left vectors */ 570 for (i=0; i<bA->nr; i++) { 571 for (j=0; j<bA->nc; j++) { 572 if (bA->m[i][j]) { 573 ierr = MatGetVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr); 574 break; 575 } 576 } 577 if (j==bA->nc) { 578 /* have an empty row */ 579 SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 580 } 581 } 582 583 ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 584 for (i=0; i<bA->nr; i++) { 585 ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 586 } 587 588 ierr = PetscFree(L);CHKERRQ(ierr); 589 } 590 PetscFunctionReturn(0); 591 } 592 593 #undef __FUNCT__ 594 #define __FUNCT__ "MatView_Nest" 595 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 596 { 597 Mat_Nest *bA = (Mat_Nest*)A->data; 598 PetscBool isascii; 599 PetscInt i,j; 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 604 if (isascii) { 605 606 PetscViewerASCIIPrintf(viewer,"Matrix object: \n"); 607 PetscViewerASCIIPushTab(viewer); /* push0 */ 608 PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 609 610 PetscViewerASCIIPrintf(viewer,"MatNest structure: \n"); 611 for (i=0; i<bA->nr; i++) { 612 for (j=0; j<bA->nc; j++) { 613 MatType type; 614 char name[256] = "",prefix[256] = ""; 615 PetscInt NR,NC; 616 PetscBool isNest = PETSC_FALSE; 617 618 if (!bA->m[i][j]) { 619 PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j); 620 continue; 621 } 622 ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 623 ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr); 624 if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 625 if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 626 ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 627 628 ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 629 630 if (isNest) { 631 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 632 ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 633 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 634 } 635 } 636 } 637 PetscViewerASCIIPopTab(viewer); /* pop0 */ 638 } 639 PetscFunctionReturn(0); 640 } 641 642 #undef __FUNCT__ 643 #define __FUNCT__ "MatZeroEntries_Nest" 644 static PetscErrorCode MatZeroEntries_Nest(Mat A) 645 { 646 Mat_Nest *bA = (Mat_Nest*)A->data; 647 PetscInt i,j; 648 PetscErrorCode ierr; 649 650 PetscFunctionBegin; 651 for (i=0; i<bA->nr; i++) { 652 for (j=0; j<bA->nc; j++) { 653 if (!bA->m[i][j]) continue; 654 ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 655 } 656 } 657 PetscFunctionReturn(0); 658 } 659 660 #undef __FUNCT__ 661 #define __FUNCT__ "MatDuplicate_Nest" 662 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 663 { 664 Mat_Nest *bA = (Mat_Nest*)A->data; 665 Mat *b; 666 PetscInt i,j,nr = bA->nr,nc = bA->nc; 667 PetscErrorCode ierr; 668 669 PetscFunctionBegin; 670 ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr); 671 for (i=0; i<nr; i++) { 672 for (j=0; j<nc; j++) { 673 if (bA->m[i][j]) { 674 ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 675 } else { 676 b[i*nc+j] = NULL; 677 } 678 } 679 } 680 ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 681 /* Give the new MatNest exclusive ownership */ 682 for (i=0; i<nr*nc; i++) { 683 ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 684 } 685 ierr = PetscFree(b);CHKERRQ(ierr); 686 687 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 688 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 689 PetscFunctionReturn(0); 690 } 691 692 /* nest api */ 693 #undef __FUNCT__ 694 #define __FUNCT__ "MatNestGetSubMat_Nest" 695 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 696 { 697 Mat_Nest *bA = (Mat_Nest*)A->data; 698 699 PetscFunctionBegin; 700 if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 701 if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 702 *mat = bA->m[idxm][jdxm]; 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "MatNestGetSubMat" 708 /*@ 709 MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 710 711 Not collective 712 713 Input Parameters: 714 + A - nest matrix 715 . idxm - index of the matrix within the nest matrix 716 - jdxm - index of the matrix within the nest matrix 717 718 Output Parameter: 719 . sub - matrix at index idxm,jdxm within the nest matrix 720 721 Level: developer 722 723 .seealso: MatNestGetSize(), MatNestGetSubMats() 724 @*/ 725 PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 726 { 727 PetscErrorCode ierr; 728 729 PetscFunctionBegin; 730 ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 731 PetscFunctionReturn(0); 732 } 733 734 #undef __FUNCT__ 735 #define __FUNCT__ "MatNestSetSubMat_Nest" 736 PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 737 { 738 Mat_Nest *bA = (Mat_Nest*)A->data; 739 PetscInt m,n,M,N,mi,ni,Mi,Ni; 740 PetscErrorCode ierr; 741 742 PetscFunctionBegin; 743 if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 744 if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 745 ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 746 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 747 ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 748 ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 749 ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 750 ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 751 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); 752 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); 753 754 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 755 ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 756 bA->m[idxm][jdxm] = mat; 757 PetscFunctionReturn(0); 758 } 759 760 #undef __FUNCT__ 761 #define __FUNCT__ "MatNestSetSubMat" 762 /*@ 763 MatNestSetSubMat - Set a single submatrix in the nest matrix. 764 765 Logically collective on the submatrix communicator 766 767 Input Parameters: 768 + A - nest matrix 769 . idxm - index of the matrix within the nest matrix 770 . jdxm - index of the matrix within the nest matrix 771 - sub - matrix at index idxm,jdxm within the nest matrix 772 773 Notes: 774 The new submatrix must have the same size and communicator as that block of the nest. 775 776 This increments the reference count of the submatrix. 777 778 Level: developer 779 780 .seealso: MatNestSetSubMats(), MatNestGetSubMat() 781 @*/ 782 PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 783 { 784 PetscErrorCode ierr; 785 786 PetscFunctionBegin; 787 ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 788 PetscFunctionReturn(0); 789 } 790 791 #undef __FUNCT__ 792 #define __FUNCT__ "MatNestGetSubMats_Nest" 793 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 794 { 795 Mat_Nest *bA = (Mat_Nest*)A->data; 796 797 PetscFunctionBegin; 798 if (M) *M = bA->nr; 799 if (N) *N = bA->nc; 800 if (mat) *mat = bA->m; 801 PetscFunctionReturn(0); 802 } 803 804 #undef __FUNCT__ 805 #define __FUNCT__ "MatNestGetSubMats" 806 /*@C 807 MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 808 809 Not collective 810 811 Input Parameters: 812 . A - nest matrix 813 814 Output Parameter: 815 + M - number of rows in the nest matrix 816 . N - number of cols in the nest matrix 817 - mat - 2d array of matrices 818 819 Notes: 820 821 The user should not free the array mat. 822 823 Level: developer 824 825 .seealso: MatNestGetSize(), MatNestGetSubMat() 826 @*/ 827 PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 828 { 829 PetscErrorCode ierr; 830 831 PetscFunctionBegin; 832 ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 833 PetscFunctionReturn(0); 834 } 835 836 #undef __FUNCT__ 837 #define __FUNCT__ "MatNestGetSize_Nest" 838 PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 839 { 840 Mat_Nest *bA = (Mat_Nest*)A->data; 841 842 PetscFunctionBegin; 843 if (M) *M = bA->nr; 844 if (N) *N = bA->nc; 845 PetscFunctionReturn(0); 846 } 847 848 #undef __FUNCT__ 849 #define __FUNCT__ "MatNestGetSize" 850 /*@ 851 MatNestGetSize - Returns the size of the nest matrix. 852 853 Not collective 854 855 Input Parameters: 856 . A - nest matrix 857 858 Output Parameter: 859 + M - number of rows in the nested mat 860 - N - number of cols in the nested mat 861 862 Notes: 863 864 Level: developer 865 866 .seealso: MatNestGetSubMat(), MatNestGetSubMats() 867 @*/ 868 PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 869 { 870 PetscErrorCode ierr; 871 872 PetscFunctionBegin; 873 ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 874 PetscFunctionReturn(0); 875 } 876 877 #undef __FUNCT__ 878 #define __FUNCT__ "MatNestGetISs_Nest" 879 static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 880 { 881 Mat_Nest *vs = (Mat_Nest*)A->data; 882 PetscInt i; 883 884 PetscFunctionBegin; 885 if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 886 if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 887 PetscFunctionReturn(0); 888 } 889 890 #undef __FUNCT__ 891 #define __FUNCT__ "MatNestGetISs" 892 /*@C 893 MatNestGetISs - Returns the index sets partitioning the row and column spaces 894 895 Not collective 896 897 Input Parameters: 898 . A - nest matrix 899 900 Output Parameter: 901 + rows - array of row index sets 902 - cols - array of column index sets 903 904 Level: advanced 905 906 Notes: 907 The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 908 909 .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs() 910 @*/ 911 PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 912 { 913 PetscErrorCode ierr; 914 915 PetscFunctionBegin; 916 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 917 ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 918 PetscFunctionReturn(0); 919 } 920 921 #undef __FUNCT__ 922 #define __FUNCT__ "MatNestGetLocalISs_Nest" 923 static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 924 { 925 Mat_Nest *vs = (Mat_Nest*)A->data; 926 PetscInt i; 927 928 PetscFunctionBegin; 929 if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 930 if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 931 PetscFunctionReturn(0); 932 } 933 934 #undef __FUNCT__ 935 #define __FUNCT__ "MatNestGetLocalISs" 936 /*@C 937 MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 938 939 Not collective 940 941 Input Parameters: 942 . A - nest matrix 943 944 Output Parameter: 945 + rows - array of row index sets (or NULL to ignore) 946 - cols - array of column index sets (or NULL to ignore) 947 948 Level: advanced 949 950 Notes: 951 The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 952 953 .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs() 954 @*/ 955 PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 956 { 957 PetscErrorCode ierr; 958 959 PetscFunctionBegin; 960 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 961 ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 962 PetscFunctionReturn(0); 963 } 964 965 #undef __FUNCT__ 966 #define __FUNCT__ "MatNestSetVecType_Nest" 967 PetscErrorCode MatNestSetVecType_Nest(Mat A,VecType vtype) 968 { 969 PetscErrorCode ierr; 970 PetscBool flg; 971 972 PetscFunctionBegin; 973 ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 974 /* In reality, this only distinguishes VECNEST and "other" */ 975 if (flg) A->ops->getvecs = MatGetVecs_Nest; 976 else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0; 977 PetscFunctionReturn(0); 978 } 979 980 #undef __FUNCT__ 981 #define __FUNCT__ "MatNestSetVecType" 982 /*@C 983 MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 984 985 Not collective 986 987 Input Parameters: 988 + A - nest matrix 989 - vtype - type to use for creating vectors 990 991 Notes: 992 993 Level: developer 994 995 .seealso: MatGetVecs() 996 @*/ 997 PetscErrorCode MatNestSetVecType(Mat A,VecType vtype) 998 { 999 PetscErrorCode ierr; 1000 1001 PetscFunctionBegin; 1002 ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr); 1003 PetscFunctionReturn(0); 1004 } 1005 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "MatNestSetSubMats_Nest" 1008 PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1009 { 1010 Mat_Nest *s = (Mat_Nest*)A->data; 1011 PetscInt i,j,m,n,M,N; 1012 PetscErrorCode ierr; 1013 1014 PetscFunctionBegin; 1015 s->nr = nr; 1016 s->nc = nc; 1017 1018 /* Create space for submatrices */ 1019 ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 1020 for (i=0; i<nr; i++) { 1021 ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 1022 } 1023 for (i=0; i<nr; i++) { 1024 for (j=0; j<nc; j++) { 1025 s->m[i][j] = a[i*nc+j]; 1026 if (a[i*nc+j]) { 1027 ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1028 } 1029 } 1030 } 1031 1032 ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1033 1034 ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 1035 ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 1036 for (i=0; i<nr; i++) s->row_len[i]=-1; 1037 for (j=0; j<nc; j++) s->col_len[j]=-1; 1038 1039 ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1040 1041 ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1042 ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1043 ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1044 ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1045 1046 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1047 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1048 1049 ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr); 1050 PetscFunctionReturn(0); 1051 } 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "MatNestSetSubMats" 1055 /*@ 1056 MatNestSetSubMats - Sets the nested submatrices 1057 1058 Collective on Mat 1059 1060 Input Parameter: 1061 + N - nested matrix 1062 . nr - number of nested row blocks 1063 . is_row - index sets for each nested row block, or NULL to make contiguous 1064 . nc - number of nested column blocks 1065 . is_col - index sets for each nested column block, or NULL to make contiguous 1066 - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1067 1068 Level: advanced 1069 1070 .seealso: MatCreateNest(), MATNEST 1071 @*/ 1072 PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1073 { 1074 PetscErrorCode ierr; 1075 PetscInt i; 1076 1077 PetscFunctionBegin; 1078 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1079 if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1080 if (nr && is_row) { 1081 PetscValidPointer(is_row,3); 1082 for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1083 } 1084 if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 1085 if (nc && is_col) { 1086 PetscValidPointer(is_col,5); 1087 for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1088 } 1089 if (nr*nc) PetscValidPointer(a,6); 1090 ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr); 1091 PetscFunctionReturn(0); 1092 } 1093 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "MatNestCreateAggregateL2G_Private" 1096 static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb) 1097 { 1098 PetscErrorCode ierr; 1099 PetscBool flg; 1100 PetscInt i,j,m,mi,*ix; 1101 1102 PetscFunctionBegin; 1103 for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 1104 if (islocal[i]) { 1105 ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 1106 flg = PETSC_TRUE; /* We found a non-trivial entry */ 1107 } else { 1108 ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 1109 } 1110 m += mi; 1111 } 1112 if (flg) { 1113 ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr); 1114 for (i=0,n=0; i<n; i++) { 1115 ISLocalToGlobalMapping smap = NULL; 1116 VecScatter scat; 1117 IS isreq; 1118 Vec lvec,gvec; 1119 union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x; 1120 Mat sub; 1121 1122 if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers"); 1123 if (colflg) { 1124 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1125 } else { 1126 ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 1127 } 1128 if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);} 1129 if (islocal[i]) { 1130 ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 1131 } else { 1132 ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 1133 } 1134 for (j=0; j<mi; j++) ix[m+j] = j; 1135 if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);} 1136 /* 1137 Now we need to extract the monolithic global indices that correspond to the given split global indices. 1138 In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces. 1139 The approach here is ugly because it uses VecScatter to move indices. 1140 */ 1141 ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr); 1142 ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr); 1143 ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr); 1144 ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr); 1145 ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 1146 for (j=0; j<mi; j++) x[j].integer = ix[m+j]; 1147 ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 1148 ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1149 ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1150 ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 1151 for (j=0; j<mi; j++) ix[m+j] = x[j].integer; 1152 ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 1153 ierr = VecDestroy(&lvec);CHKERRQ(ierr); 1154 ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1155 ierr = ISDestroy(&isreq);CHKERRQ(ierr); 1156 ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 1157 m += mi; 1158 } 1159 ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 1160 *ltogb = NULL; 1161 } else { 1162 *ltog = NULL; 1163 *ltogb = NULL; 1164 } 1165 PetscFunctionReturn(0); 1166 } 1167 1168 1169 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1170 /* 1171 nprocessors = NP 1172 Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1)) 1173 proc 0: => (g_0,h_0,) 1174 proc 1: => (g_1,h_1,) 1175 ... 1176 proc nprocs-1: => (g_NP-1,h_NP-1,) 1177 1178 proc 0: proc 1: proc nprocs-1: 1179 is[0] = (0,1,2,...,nlocal(g_0)-1) (0,1,...,nlocal(g_1)-1) (0,1,...,nlocal(g_NP-1)) 1180 1181 proc 0: 1182 is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1) 1183 proc 1: 1184 is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1) 1185 1186 proc NP-1: 1187 is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1) 1188 */ 1189 #undef __FUNCT__ 1190 #define __FUNCT__ "MatSetUp_NestIS_Private" 1191 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1192 { 1193 Mat_Nest *vs = (Mat_Nest*)A->data; 1194 PetscInt i,j,offset,n,nsum,bs; 1195 PetscErrorCode ierr; 1196 Mat sub = NULL; 1197 1198 PetscFunctionBegin; 1199 ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 1200 ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1201 if (is_row) { /* valid IS is passed in */ 1202 /* refs on is[] are incremeneted */ 1203 for (i=0; i<vs->nr; i++) { 1204 ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1205 1206 vs->isglobal.row[i] = is_row[i]; 1207 } 1208 } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 1209 nsum = 0; 1210 for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 1211 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1212 if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 1213 ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 1214 if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 1215 nsum += n; 1216 } 1217 ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 1218 offset -= nsum; 1219 for (i=0; i<vs->nr; i++) { 1220 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1221 ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 1222 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1223 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1224 ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 1225 offset += n; 1226 } 1227 } 1228 1229 if (is_col) { /* valid IS is passed in */ 1230 /* refs on is[] are incremeneted */ 1231 for (j=0; j<vs->nc; j++) { 1232 ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1233 1234 vs->isglobal.col[j] = is_col[j]; 1235 } 1236 } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 1237 offset = A->cmap->rstart; 1238 nsum = 0; 1239 for (j=0; j<vs->nc; j++) { 1240 ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1241 if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 1242 ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 1243 if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 1244 nsum += n; 1245 } 1246 ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 1247 offset -= nsum; 1248 for (j=0; j<vs->nc; j++) { 1249 ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1250 ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 1251 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1252 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1253 ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 1254 offset += n; 1255 } 1256 } 1257 1258 /* Set up the local ISs */ 1259 ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 1260 ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 1261 for (i=0,offset=0; i<vs->nr; i++) { 1262 IS isloc; 1263 ISLocalToGlobalMapping rmap = NULL; 1264 PetscInt nlocal,bs; 1265 ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1266 if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);} 1267 if (rmap) { 1268 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1269 ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1270 ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1271 ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1272 } else { 1273 nlocal = 0; 1274 isloc = NULL; 1275 } 1276 vs->islocal.row[i] = isloc; 1277 offset += nlocal; 1278 } 1279 for (i=0,offset=0; i<vs->nc; i++) { 1280 IS isloc; 1281 ISLocalToGlobalMapping cmap = NULL; 1282 PetscInt nlocal,bs; 1283 ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 1284 if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);} 1285 if (cmap) { 1286 ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1287 ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1288 ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1289 ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1290 } else { 1291 nlocal = 0; 1292 isloc = NULL; 1293 } 1294 vs->islocal.col[i] = isloc; 1295 offset += nlocal; 1296 } 1297 1298 /* Set up the aggregate ISLocalToGlobalMapping */ 1299 { 1300 ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb; 1301 ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr); 1302 ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr); 1303 if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 1304 if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);} 1305 ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 1306 ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr); 1307 ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 1308 ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr); 1309 } 1310 1311 #if defined(PETSC_USE_DEBUG) 1312 for (i=0; i<vs->nr; i++) { 1313 for (j=0; j<vs->nc; j++) { 1314 PetscInt m,n,M,N,mi,ni,Mi,Ni; 1315 Mat B = vs->m[i][j]; 1316 if (!B) continue; 1317 ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 1318 ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 1319 ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 1320 ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 1321 ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 1322 ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 1323 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); 1324 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); 1325 } 1326 } 1327 #endif 1328 1329 /* Set A->assembled if all non-null blocks are currently assembled */ 1330 for (i=0; i<vs->nr; i++) { 1331 for (j=0; j<vs->nc; j++) { 1332 if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1333 } 1334 } 1335 A->assembled = PETSC_TRUE; 1336 PetscFunctionReturn(0); 1337 } 1338 1339 #undef __FUNCT__ 1340 #define __FUNCT__ "MatCreateNest" 1341 /*@C 1342 MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1343 1344 Collective on Mat 1345 1346 Input Parameter: 1347 + comm - Communicator for the new Mat 1348 . nr - number of nested row blocks 1349 . is_row - index sets for each nested row block, or NULL to make contiguous 1350 . nc - number of nested column blocks 1351 . is_col - index sets for each nested column block, or NULL to make contiguous 1352 - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1353 1354 Output Parameter: 1355 . B - new matrix 1356 1357 Level: advanced 1358 1359 .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST 1360 @*/ 1361 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1362 { 1363 Mat A; 1364 PetscErrorCode ierr; 1365 1366 PetscFunctionBegin; 1367 *B = 0; 1368 ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1369 ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 1370 ierr = MatSetUp(A);CHKERRQ(ierr); 1371 ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1372 *B = A; 1373 PetscFunctionReturn(0); 1374 } 1375 1376 #undef __FUNCT__ 1377 #define __FUNCT__ "MatConvert_Nest_AIJ" 1378 PETSC_EXTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1379 { 1380 PetscErrorCode ierr; 1381 Mat_Nest *nest = (Mat_Nest*)A->data; 1382 PetscInt m,n,M,N,i,j,k,*dnnz,*onnz; 1383 Mat C; 1384 1385 PetscFunctionBegin; 1386 ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1387 ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1388 switch (reuse) { 1389 case MAT_INITIAL_MATRIX: 1390 ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 1391 ierr = MatSetType(C,newtype);CHKERRQ(ierr); 1392 ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 1393 *newmat = C; 1394 break; 1395 case MAT_REUSE_MATRIX: 1396 C = *newmat; 1397 break; 1398 default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse"); 1399 } 1400 1401 /* Preallocation */ 1402 ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr); 1403 onnz = dnnz + m; 1404 for (k=0; k<m; k++) { 1405 dnnz[k] = 0; 1406 onnz[k] = 0; 1407 } 1408 for (j=0; j<nest->nc; ++j) { 1409 IS bNis; 1410 PetscInt bN; 1411 const PetscInt *bNindices; 1412 /* Using global column indices and ISAllGather() is not scalable. */ 1413 ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1414 ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr); 1415 ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1416 for (i=0; i<nest->nr; ++i) { 1417 PetscSF bmsf; 1418 PetscSFNode *bmedges; 1419 Mat B; 1420 PetscInt bm, *bmdnnz, br; 1421 const PetscInt *bmindices; 1422 B = nest->m[i][j]; 1423 if (!B) continue; 1424 ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1425 ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1426 ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr); 1427 ierr = PetscMalloc1(2*bm,&bmedges);CHKERRQ(ierr); 1428 ierr = PetscMalloc1(2*bm,&bmdnnz);CHKERRQ(ierr); 1429 for (k = 0; k < 2*bm; ++k) bmdnnz[k] = 0; 1430 /* 1431 Locate the owners for all of the locally-owned global row indices for this row block. 1432 These determine the roots of PetscSF used to communicate preallocation data to row owners. 1433 The roots correspond to the dnnz and onnz entries; thus, there are two roots per row. 1434 */ 1435 for (br = 0; br < bm; ++br) { 1436 PetscInt row = bmindices[br], rowowner = 0, brncols, col, colowner = 0; 1437 const PetscInt *brcols; 1438 PetscInt rowrel = 0; /* row's relative index on its owner rank */ 1439 PetscInt rowownerm; /* local row size on row's owning rank. */ 1440 1441 ierr = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr); 1442 rowownerm = A->rmap->range[rowowner+1]-A->rmap->range[rowowner]; 1443 1444 bmedges[br].rank = rowowner; bmedges[br].index = rowrel; /* edge from bmdnnz to dnnz */ 1445 bmedges[br].rank = rowowner; bmedges[br].index = rowrel+rowownerm; /* edge from bmonnz to onnz */ 1446 /* Now actually compute the data -- bmdnnz and bmonnz by looking at the global columns in the br row of this block. */ 1447 /* Note that this is not a pessimistic bound only because we assume the index sets embedding the blocks do not overlap. */ 1448 ierr = MatGetRow(B,br,&brncols,&brcols,NULL);CHKERRQ(ierr); 1449 for (k=0; k<brncols; k++) { 1450 col = bNindices[brcols[k]]; 1451 ierr = PetscLayoutFindOwnerIndex(A->cmap,col,&colowner,NULL);CHKERRQ(ierr); 1452 if (colowner == rowowner) bmdnnz[br]++; 1453 else onnz[br]++; 1454 } 1455 ierr = MatRestoreRow(B,br,&brncols,&brcols,NULL);CHKERRQ(ierr); 1456 } 1457 ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1458 /* bsf will have to take care of disposing of bedges. */ 1459 ierr = PetscSFSetGraph(bmsf,m,2*bm,NULL,PETSC_COPY_VALUES,bmedges,PETSC_OWN_POINTER);CHKERRQ(ierr); 1460 ierr = PetscSFReduceBegin(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr); 1461 ierr = PetscSFReduceEnd(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr); 1462 ierr = PetscFree(bmdnnz);CHKERRQ(ierr); 1463 ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr); 1464 } 1465 ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1466 ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1467 } 1468 ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr); 1469 ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr); 1470 ierr = PetscFree(dnnz);CHKERRQ(ierr); 1471 1472 /* Fill by row */ 1473 for (j=0; j<nest->nc; ++j) { 1474 /* Using global column indices and ISAllGather() is not scalable. */ 1475 IS bNis; 1476 PetscInt bN; 1477 const PetscInt *bNindices; 1478 ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1479 ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr); 1480 ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1481 for (i=0; i<nest->nr; ++i) { 1482 Mat B; 1483 PetscInt bm, br; 1484 const PetscInt *bmindices; 1485 B = nest->m[i][j]; 1486 if (!B) continue; 1487 ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1488 ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1489 for (br = 0; br < bm; ++br) { 1490 PetscInt row = bmindices[br], brncols, *cols; 1491 const PetscInt *brcols; 1492 const PetscScalar *brcoldata; 1493 ierr = MatGetRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1494 ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr); 1495 for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]]; 1496 /* 1497 Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match. 1498 Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES. 1499 */ 1500 ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr); 1501 ierr = MatRestoreRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1502 ierr = PetscFree(cols);CHKERRQ(ierr); 1503 } 1504 ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1505 } 1506 ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1507 ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1508 } 1509 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1510 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1511 PetscFunctionReturn(0); 1512 } 1513 1514 /*MC 1515 MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1516 1517 Level: intermediate 1518 1519 Notes: 1520 This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1521 It allows the use of symmetric and block formats for parts of multi-physics simulations. 1522 It is usually used with DMComposite and DMCreateMatrix() 1523 1524 .seealso: MatCreate(), MatType, MatCreateNest() 1525 M*/ 1526 #undef __FUNCT__ 1527 #define __FUNCT__ "MatCreate_Nest" 1528 PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A) 1529 { 1530 Mat_Nest *s; 1531 PetscErrorCode ierr; 1532 1533 PetscFunctionBegin; 1534 ierr = PetscNewLog(A,&s);CHKERRQ(ierr); 1535 A->data = (void*)s; 1536 1537 s->nr = -1; 1538 s->nc = -1; 1539 s->m = NULL; 1540 s->splitassembly = PETSC_FALSE; 1541 1542 ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 1543 1544 A->ops->mult = MatMult_Nest; 1545 A->ops->multadd = MatMultAdd_Nest; 1546 A->ops->multtranspose = MatMultTranspose_Nest; 1547 A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1548 A->ops->assemblybegin = MatAssemblyBegin_Nest; 1549 A->ops->assemblyend = MatAssemblyEnd_Nest; 1550 A->ops->zeroentries = MatZeroEntries_Nest; 1551 A->ops->duplicate = MatDuplicate_Nest; 1552 A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1553 A->ops->destroy = MatDestroy_Nest; 1554 A->ops->view = MatView_Nest; 1555 A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1556 A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1557 A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1558 A->ops->getdiagonal = MatGetDiagonal_Nest; 1559 A->ops->diagonalscale = MatDiagonalScale_Nest; 1560 A->ops->scale = MatScale_Nest; 1561 A->ops->shift = MatShift_Nest; 1562 1563 A->spptr = 0; 1564 A->same_nonzero = PETSC_FALSE; 1565 A->assembled = PETSC_FALSE; 1566 1567 /* expose Nest api's */ 1568 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C", MatNestGetSubMat_Nest);CHKERRQ(ierr); 1569 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1570 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1571 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C", MatNestGetSize_Nest);CHKERRQ(ierr); 1572 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C", MatNestGetISs_Nest);CHKERRQ(ierr); 1573 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 1574 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C", MatNestSetVecType_Nest);CHKERRQ(ierr); 1575 ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1576 1577 ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1578 PetscFunctionReturn(0); 1579 } 1580