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