1d8588912SDave May 2ec9811eeSJed Brown #include "matnestimpl.h" /*I "petscmat.h" I*/ 3d8588912SDave May 4c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]); 57874fa86SDave May static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left); 6c8883902SJed Brown 7d8588912SDave May /* private functions */ 8d8588912SDave May #undef __FUNCT__ 98188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private" 108188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N) 11d8588912SDave May { 12d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 138188e55aSJed Brown PetscInt i,j; 14d8588912SDave May PetscErrorCode ierr; 15d8588912SDave May 16d8588912SDave May PetscFunctionBegin; 178188e55aSJed Brown *m = *n = *M = *N = 0; 188188e55aSJed Brown for (i=0; i<bA->nr; i++) { /* rows */ 198188e55aSJed Brown PetscInt sm,sM; 208188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr); 218188e55aSJed Brown ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr); 228188e55aSJed Brown *m += sm; 238188e55aSJed Brown *M += sM; 24d8588912SDave May } 258188e55aSJed Brown for (j=0; j<bA->nc; j++) { /* cols */ 268188e55aSJed Brown PetscInt sn,sN; 278188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr); 288188e55aSJed Brown ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr); 298188e55aSJed Brown *n += sn; 308188e55aSJed Brown *N += sN; 31d8588912SDave May } 32d8588912SDave May PetscFunctionReturn(0); 33d8588912SDave May } 34d8588912SDave May 35d8588912SDave May #undef __FUNCT__ 36d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2" 37207556f9SJed Brown PETSC_UNUSED 38d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y) 39d8588912SDave May { 40d8588912SDave May PetscBool isANest, isxNest, isyNest; 41d8588912SDave May PetscErrorCode ierr; 42d8588912SDave May 43d8588912SDave May PetscFunctionBegin; 44d8588912SDave May isANest = isxNest = isyNest = PETSC_FALSE; 45d8588912SDave May ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr); 46d8588912SDave May ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr); 47d8588912SDave May ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr); 48d8588912SDave May 49d8588912SDave May if (isANest && isxNest && isyNest) { 50d8588912SDave May PetscFunctionReturn(0); 51d8588912SDave May } else { 52d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations"); 53d8588912SDave May PetscFunctionReturn(0); 54d8588912SDave May } 55d8588912SDave May PetscFunctionReturn(0); 56d8588912SDave May } 57d8588912SDave May 58d8588912SDave May /* 59d8588912SDave May This function is called every time we insert a sub matrix the Nest. 60d8588912SDave May It ensures that every Nest along row r, and col c has the same dimensions 61d8588912SDave May as the submat being inserted. 62d8588912SDave May */ 63d8588912SDave May #undef __FUNCT__ 64d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckConsistency" 65063a28d4SJed Brown PETSC_UNUSED 66d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c) 67d8588912SDave May { 68d8588912SDave May Mat_Nest *b = (Mat_Nest*)A->data; 69d8588912SDave May PetscInt i,j; 70d8588912SDave May PetscInt nr,nc; 71d8588912SDave May PetscInt sM,sN,M,N; 72d8588912SDave May Mat mat; 73d8588912SDave May PetscErrorCode ierr; 74d8588912SDave May 75d8588912SDave May PetscFunctionBegin; 76d8588912SDave May nr = b->nr; 77d8588912SDave May nc = b->nc; 78d8588912SDave May ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr); 79d8588912SDave May for (i=0; i<nr; i++) { 80d8588912SDave May mat = b->m[i][c]; 81d8588912SDave May if (mat) { 82d8588912SDave May ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 83d8588912SDave May /* Check columns match */ 84d8588912SDave May if (sN != N) { 85d8588912SDave May SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N ); 86d8588912SDave May } 87d8588912SDave May } 88d8588912SDave May } 89d8588912SDave May 90d8588912SDave May for (j=0; j<nc; j++) { 91d8588912SDave May mat = b->m[r][j]; 92d8588912SDave May if (mat) { 93d8588912SDave May ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 94d8588912SDave May /* Check rows match */ 95d8588912SDave May if (sM != M) { 96d8588912SDave May SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M ); 97d8588912SDave May } 98d8588912SDave May } 99d8588912SDave May } 100d8588912SDave May PetscFunctionReturn(0); 101d8588912SDave May } 102d8588912SDave May 103d8588912SDave May /* 104d8588912SDave May Checks if the row/col's contain a complete set of IS's. 105d8588912SDave May Once they are filled, the offsets are computed. This saves having to update 106d8588912SDave May the index set entries each time we insert something new. 107d8588912SDave May */ 108d8588912SDave May #undef __FUNCT__ 109d8588912SDave May #define __FUNCT__ "PETSc_MatNest_UpdateStructure" 110063a28d4SJed Brown PETSC_UNUSED 111d8588912SDave May static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx) 112d8588912SDave May { 113d8588912SDave May Mat_Nest *b = (Mat_Nest*)A->data; 114d8588912SDave May PetscInt m,n,i,j; 115d8588912SDave May PetscBool fullrow,fullcol; 116d8588912SDave May PetscErrorCode ierr; 117d8588912SDave May 118d8588912SDave May PetscFunctionBegin; 119d8588912SDave May ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr); 120d8588912SDave May b->row_len[ridx] = m; 121d8588912SDave May b->col_len[cidx] = n; 122d8588912SDave May 123d8588912SDave May fullrow = PETSC_TRUE; 124d8588912SDave May for (i=0; i<b->nr; i++) { 125d8588912SDave May if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; } 126d8588912SDave May } 127f349c1fdSJed Brown if ( (fullrow) && (!b->isglobal.row[0]) ){ 128d8588912SDave May PetscInt cnt = 0; 129d8588912SDave May for (i=0; i<b->nr; i++) { 130f349c1fdSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr); 131d8588912SDave May cnt = cnt + b->row_len[i]; 132d8588912SDave May } 133d8588912SDave May } 134d8588912SDave May 135d8588912SDave May fullcol = PETSC_TRUE; 136d8588912SDave May for (j=0; j<b->nc; j++) { 137d8588912SDave May if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; } 138d8588912SDave May } 139f349c1fdSJed Brown if( (fullcol) && (!b->isglobal.col[0]) ){ 140d8588912SDave May PetscInt cnt = 0; 141d8588912SDave May for (j=0; j<b->nc; j++) { 142f349c1fdSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr); 143d8588912SDave May cnt = cnt + b->col_len[j]; 144d8588912SDave May } 145d8588912SDave May } 146d8588912SDave May PetscFunctionReturn(0); 147d8588912SDave May } 148d8588912SDave May 149d8588912SDave May /* operations */ 150d8588912SDave May #undef __FUNCT__ 151d8588912SDave May #define __FUNCT__ "MatMult_Nest" 152207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 153d8588912SDave May { 154d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 155207556f9SJed Brown Vec *bx = bA->right,*by = bA->left; 156207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 157d8588912SDave May PetscErrorCode ierr; 158d8588912SDave May 159d8588912SDave May PetscFunctionBegin; 160207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 161207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 162207556f9SJed Brown for (i=0; i<nr; i++) { 163d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 164207556f9SJed Brown for (j=0; j<nc; j++) { 165207556f9SJed Brown if (!bA->m[i][j]) continue; 166d8588912SDave May /* y[i] <- y[i] + A[i][j] * x[j] */ 167d8588912SDave May ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 168d8588912SDave May } 169d8588912SDave May } 170207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 171207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 172d8588912SDave May PetscFunctionReturn(0); 173d8588912SDave May } 174d8588912SDave May 175d8588912SDave May #undef __FUNCT__ 176d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest" 177207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 178d8588912SDave May { 179d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 180207556f9SJed Brown Vec *bx = bA->left,*by = bA->right; 181207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 182d8588912SDave May PetscErrorCode ierr; 183d8588912SDave May 184d8588912SDave May PetscFunctionBegin; 185d8588912SDave May if (A->symmetric) { 186d8588912SDave May ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr); 187d8588912SDave May PetscFunctionReturn(0); 188d8588912SDave May } 189207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 190207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 191207556f9SJed Brown for (i=0; i<nr; i++) { 192d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 193207556f9SJed Brown for (j=0; j<nc; j++) { 194207556f9SJed Brown if (!bA->m[j][i]) continue; 195d8588912SDave May /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */ 196d8588912SDave May ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr); 197d8588912SDave May } 198d8588912SDave May } 199207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 200207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 201d8588912SDave May PetscFunctionReturn(0); 202d8588912SDave May } 203d8588912SDave May 204d8588912SDave May #undef __FUNCT__ 205e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList" 206e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 207e2d7f03fSJed Brown { 208e2d7f03fSJed Brown PetscErrorCode ierr; 209e2d7f03fSJed Brown IS *lst = *list; 210e2d7f03fSJed Brown PetscInt i; 211e2d7f03fSJed Brown 212e2d7f03fSJed Brown PetscFunctionBegin; 213e2d7f03fSJed Brown if (!lst) PetscFunctionReturn(0); 2146bf464f9SBarry Smith for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);} 215e2d7f03fSJed Brown ierr = PetscFree(lst);CHKERRQ(ierr); 216e2d7f03fSJed Brown *list = PETSC_NULL; 217e2d7f03fSJed Brown PetscFunctionReturn(0); 218e2d7f03fSJed Brown } 219e2d7f03fSJed Brown 220e2d7f03fSJed Brown #undef __FUNCT__ 221d8588912SDave May #define __FUNCT__ "MatDestroy_Nest" 222207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A) 223d8588912SDave May { 224d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 225d8588912SDave May PetscInt i,j; 226d8588912SDave May PetscErrorCode ierr; 227d8588912SDave May 228d8588912SDave May PetscFunctionBegin; 229d8588912SDave May /* release the matrices and the place holders */ 230e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 231e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 232e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 233e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 234d8588912SDave May 235d8588912SDave May ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 236d8588912SDave May ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 237d8588912SDave May 238207556f9SJed Brown ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 239207556f9SJed Brown 240d8588912SDave May /* release the matrices and the place holders */ 241d8588912SDave May if (vs->m) { 242d8588912SDave May for (i=0; i<vs->nr; i++) { 243d8588912SDave May for (j=0; j<vs->nc; j++) { 2446bf464f9SBarry Smith ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr); 245d8588912SDave May } 246d8588912SDave May ierr = PetscFree( vs->m[i] );CHKERRQ(ierr); 247d8588912SDave May } 248d8588912SDave May ierr = PetscFree(vs->m);CHKERRQ(ierr); 249d8588912SDave May } 250bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 251d8588912SDave May 252c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "",0);CHKERRQ(ierr); 253*0782ca92SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "",0);CHKERRQ(ierr); 254c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "",0);CHKERRQ(ierr); 255c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "",0);CHKERRQ(ierr); 256c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "",0);CHKERRQ(ierr); 257c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "",0);CHKERRQ(ierr); 258d8588912SDave May PetscFunctionReturn(0); 259d8588912SDave May } 260d8588912SDave May 261d8588912SDave May #undef __FUNCT__ 262d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest" 263207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 264d8588912SDave May { 265d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 266d8588912SDave May PetscInt i,j; 267d8588912SDave May PetscErrorCode ierr; 268d8588912SDave May 269d8588912SDave May PetscFunctionBegin; 270d8588912SDave May for (i=0; i<vs->nr; i++) { 271d8588912SDave May for (j=0; j<vs->nc; j++) { 272d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); } 273d8588912SDave May } 274d8588912SDave May } 275d8588912SDave May PetscFunctionReturn(0); 276d8588912SDave May } 277d8588912SDave May 278d8588912SDave May #undef __FUNCT__ 279d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 280207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 281d8588912SDave May { 282d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 283d8588912SDave May PetscInt i,j; 284d8588912SDave May PetscErrorCode ierr; 285d8588912SDave May 286d8588912SDave May PetscFunctionBegin; 287d8588912SDave May for (i=0; i<vs->nr; i++) { 288d8588912SDave May for (j=0; j<vs->nc; j++) { 289d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); } 290d8588912SDave May } 291d8588912SDave May } 292d8588912SDave May PetscFunctionReturn(0); 293d8588912SDave May } 294d8588912SDave May 295d8588912SDave May #undef __FUNCT__ 296f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 297f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 298d8588912SDave May { 299207556f9SJed Brown PetscErrorCode ierr; 300f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 301f349c1fdSJed Brown PetscInt j; 302f349c1fdSJed Brown Mat sub; 303d8588912SDave May 304d8588912SDave May PetscFunctionBegin; 3058188e55aSJed Brown sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */ 306f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 3078188e55aSJed Brown if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 308f349c1fdSJed Brown *B = sub; 309f349c1fdSJed Brown PetscFunctionReturn(0); 310d8588912SDave May } 311d8588912SDave May 312f349c1fdSJed Brown #undef __FUNCT__ 313f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 314f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 315f349c1fdSJed Brown { 316207556f9SJed Brown PetscErrorCode ierr; 317f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 318f349c1fdSJed Brown PetscInt i; 319f349c1fdSJed Brown Mat sub; 320f349c1fdSJed Brown 321f349c1fdSJed Brown PetscFunctionBegin; 3228188e55aSJed Brown sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */ 323f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 3248188e55aSJed Brown if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 325f349c1fdSJed Brown *B = sub; 326f349c1fdSJed Brown PetscFunctionReturn(0); 327d8588912SDave May } 328d8588912SDave May 329f349c1fdSJed Brown #undef __FUNCT__ 330f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 331f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 332f349c1fdSJed Brown { 333f349c1fdSJed Brown PetscErrorCode ierr; 334f349c1fdSJed Brown PetscInt i; 335f349c1fdSJed Brown PetscBool flg; 336f349c1fdSJed Brown 337f349c1fdSJed Brown PetscFunctionBegin; 338f349c1fdSJed Brown PetscValidPointer(list,3); 339f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 340f349c1fdSJed Brown PetscValidIntPointer(found,5); 341f349c1fdSJed Brown *found = -1; 342f349c1fdSJed Brown for (i=0; i<n; i++) { 343207556f9SJed Brown if (!list[i]) continue; 344f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 345f349c1fdSJed Brown if (flg) { 346f349c1fdSJed Brown *found = i; 347f349c1fdSJed Brown PetscFunctionReturn(0); 348f349c1fdSJed Brown } 349f349c1fdSJed Brown } 350f349c1fdSJed Brown SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set"); 351f349c1fdSJed Brown PetscFunctionReturn(0); 352f349c1fdSJed Brown } 353f349c1fdSJed Brown 354f349c1fdSJed Brown #undef __FUNCT__ 3558188e55aSJed Brown #define __FUNCT__ "MatNestGetRow" 3568188e55aSJed Brown /* Get a block row as a new MatNest */ 3578188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 3588188e55aSJed Brown { 3598188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3608188e55aSJed Brown Mat C; 3618188e55aSJed Brown char keyname[256]; 3628188e55aSJed Brown PetscErrorCode ierr; 3638188e55aSJed Brown 3648188e55aSJed Brown PetscFunctionBegin; 3658188e55aSJed Brown *B = PETSC_NULL; 3668188e55aSJed Brown ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr); 3678188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 3688188e55aSJed Brown if (*B) PetscFunctionReturn(0); 3698188e55aSJed Brown 3708188e55aSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 3718188e55aSJed Brown (*B)->assembled = A->assembled; 3728188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 3738188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 3748188e55aSJed Brown PetscFunctionReturn(0); 3758188e55aSJed Brown } 3768188e55aSJed Brown 3778188e55aSJed Brown #undef __FUNCT__ 378f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 379f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 380f349c1fdSJed Brown { 381f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3828188e55aSJed Brown PetscErrorCode ierr; 3838188e55aSJed Brown PetscInt row,col,i; 3848188e55aSJed Brown IS *iscopy; 3858188e55aSJed Brown Mat *matcopy; 3868188e55aSJed Brown PetscBool same,isFullCol; 387f349c1fdSJed Brown 388f349c1fdSJed Brown PetscFunctionBegin; 3898188e55aSJed Brown /* Check if full column space. This is a hack */ 3908188e55aSJed Brown isFullCol = PETSC_FALSE; 3918188e55aSJed Brown ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 3928188e55aSJed Brown if (same) { 3938188e55aSJed Brown PetscInt n,first,step; 3948188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 3958188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 3968188e55aSJed Brown if (A->cmap->n == n && A->cmap->rstart == first && step == 1) isFullCol = PETSC_TRUE; 3978188e55aSJed Brown } 3988188e55aSJed Brown 3998188e55aSJed Brown if (isFullCol) { 4008188e55aSJed Brown PetscInt row; 4018188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 4028188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 4038188e55aSJed Brown } else { 404f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 405f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 406f349c1fdSJed Brown *B = vs->m[row][col]; 4078188e55aSJed Brown } 408f349c1fdSJed Brown PetscFunctionReturn(0); 409f349c1fdSJed Brown } 410f349c1fdSJed Brown 411f349c1fdSJed Brown #undef __FUNCT__ 412f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 413f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 414f349c1fdSJed Brown { 415f349c1fdSJed Brown PetscErrorCode ierr; 416f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 417f349c1fdSJed Brown Mat sub; 418f349c1fdSJed Brown 419f349c1fdSJed Brown PetscFunctionBegin; 420f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 421f349c1fdSJed Brown switch (reuse) { 422f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 4237874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 424f349c1fdSJed Brown *B = sub; 425f349c1fdSJed Brown break; 426f349c1fdSJed Brown case MAT_REUSE_MATRIX: 427f349c1fdSJed Brown if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 428f349c1fdSJed Brown break; 429f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 430f349c1fdSJed Brown break; 431f349c1fdSJed Brown } 432f349c1fdSJed Brown PetscFunctionReturn(0); 433f349c1fdSJed Brown } 434f349c1fdSJed Brown 435f349c1fdSJed Brown #undef __FUNCT__ 436f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 437f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 438f349c1fdSJed Brown { 439f349c1fdSJed Brown PetscErrorCode ierr; 440f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 441f349c1fdSJed Brown Mat sub; 442f349c1fdSJed Brown 443f349c1fdSJed Brown PetscFunctionBegin; 444f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 445f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 446f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 447f349c1fdSJed Brown *B = sub; 448d8588912SDave May PetscFunctionReturn(0); 449d8588912SDave May } 450d8588912SDave May 451d8588912SDave May #undef __FUNCT__ 452d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 453207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 454d8588912SDave May { 455d8588912SDave May PetscErrorCode ierr; 456f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 457f349c1fdSJed Brown Mat sub; 458d8588912SDave May 459d8588912SDave May PetscFunctionBegin; 460f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 461f349c1fdSJed Brown if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 462f349c1fdSJed Brown if (sub) { 463f349c1fdSJed Brown if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4646bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 465d8588912SDave May } 466d8588912SDave May PetscFunctionReturn(0); 467d8588912SDave May } 468d8588912SDave May 469d8588912SDave May #undef __FUNCT__ 4707874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 4717874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 4727874fa86SDave May { 4737874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4747874fa86SDave May Vec *bdiag; 4757874fa86SDave May PetscInt i; 4767874fa86SDave May PetscErrorCode ierr; 4777874fa86SDave May 4787874fa86SDave May PetscFunctionBegin; 479a0769a71SSatish Balay /* ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); */ 480a0769a71SSatish Balay /* ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); */ 4817874fa86SDave May ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr); 4827874fa86SDave May for (i=0; i<bA->nr; i++) { 4837874fa86SDave May if (bA->m[i][i]) { 4847874fa86SDave May ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr); 4857874fa86SDave May } else { 4867874fa86SDave May ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr); 4877874fa86SDave May } 4887874fa86SDave May } 4897874fa86SDave May PetscFunctionReturn(0); 4907874fa86SDave May } 4917874fa86SDave May 4927874fa86SDave May #undef __FUNCT__ 4937874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest2" 4947874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v) 4957874fa86SDave May { 4967874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4977874fa86SDave May Vec diag,*bdiag; 4987874fa86SDave May VecScatter *vscat; 4997874fa86SDave May PetscInt i; 5007874fa86SDave May PetscErrorCode ierr; 5017874fa86SDave May 5027874fa86SDave May PetscFunctionBegin; 5037874fa86SDave May ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); 5047874fa86SDave May ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr); 5057874fa86SDave May 5067874fa86SDave May /* scatter diag into v */ 5077874fa86SDave May ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); 5087874fa86SDave May ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr); 5097874fa86SDave May for (i=0; i<bA->nr; i++) { 5107874fa86SDave May ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr); 5117874fa86SDave May ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5127874fa86SDave May } 5137874fa86SDave May for (i=0; i<bA->nr; i++) { 5147874fa86SDave May ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5157874fa86SDave May } 5167874fa86SDave May for (i=0; i<bA->nr; i++) { 5176bf464f9SBarry Smith ierr = VecScatterDestroy(&vscat[i]);CHKERRQ(ierr); 5187874fa86SDave May } 5197874fa86SDave May ierr = PetscFree(vscat);CHKERRQ(ierr); 5206bf464f9SBarry Smith ierr = VecDestroy(&diag);CHKERRQ(ierr); 5217874fa86SDave May PetscFunctionReturn(0); 5227874fa86SDave May } 5237874fa86SDave May 5247874fa86SDave May #undef __FUNCT__ 5257874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 5267874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 5277874fa86SDave May { 5287874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 5297874fa86SDave May Vec *bl,*br; 5307874fa86SDave May PetscInt i,j; 5317874fa86SDave May PetscErrorCode ierr; 5327874fa86SDave May 5337874fa86SDave May PetscFunctionBegin; 5347874fa86SDave May ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr); 5357874fa86SDave May ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr); 5367874fa86SDave May for (i=0; i<bA->nr; i++) { 5377874fa86SDave May for (j=0; j<bA->nc; j++) { 5387874fa86SDave May if (bA->m[i][j]) { 5397874fa86SDave May ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr); 5407874fa86SDave May } 5417874fa86SDave May } 5427874fa86SDave May } 5437874fa86SDave May PetscFunctionReturn(0); 5447874fa86SDave May } 5457874fa86SDave May 5467874fa86SDave May #undef __FUNCT__ 5477874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest2" 5487874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r) 5497874fa86SDave May { 5507874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 5517874fa86SDave May Vec bl,br,*ble,*bre; 5527874fa86SDave May VecScatter *vscatl,*vscatr; 5537874fa86SDave May PetscInt i; 5547874fa86SDave May PetscErrorCode ierr; 5557874fa86SDave May 5567874fa86SDave May PetscFunctionBegin; 5577874fa86SDave May /* scatter l,r into bl,br */ 5587874fa86SDave May ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr); 5597874fa86SDave May ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr); 5607874fa86SDave May ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr); 5617874fa86SDave May 5627874fa86SDave May /* row */ 5637874fa86SDave May ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr); 5647874fa86SDave May for (i=0; i<bA->nr; i++) { 5657874fa86SDave May ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr); 5667874fa86SDave May ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5677874fa86SDave May } 5687874fa86SDave May for (i=0; i<bA->nr; i++) { 5697874fa86SDave May ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5707874fa86SDave May } 5717874fa86SDave May /* col */ 5727874fa86SDave May ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr); 5737874fa86SDave May for (i=0; i<bA->nc; i++) { 5747874fa86SDave May ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr); 5757874fa86SDave May ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5767874fa86SDave May } 5777874fa86SDave May for (i=0; i<bA->nc; i++) { 5787874fa86SDave May ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5797874fa86SDave May } 5807874fa86SDave May 5817874fa86SDave May ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr); 5827874fa86SDave May 5837874fa86SDave May for (i=0; i<bA->nr; i++) { 5846bf464f9SBarry Smith ierr = VecScatterDestroy(&vscatl[i]);CHKERRQ(ierr); 5857874fa86SDave May } 5867874fa86SDave May for (i=0; i<bA->nc; i++) { 5876bf464f9SBarry Smith ierr = VecScatterDestroy(&vscatr[i]);CHKERRQ(ierr); 5887874fa86SDave May } 5897874fa86SDave May ierr = PetscFree(vscatl);CHKERRQ(ierr); 5907874fa86SDave May ierr = PetscFree(vscatr);CHKERRQ(ierr); 5916bf464f9SBarry Smith ierr = VecDestroy(&bl);CHKERRQ(ierr); 5926bf464f9SBarry Smith ierr = VecDestroy(&br);CHKERRQ(ierr); 5937874fa86SDave May PetscFunctionReturn(0); 5947874fa86SDave May } 5957874fa86SDave May 5967874fa86SDave May #undef __FUNCT__ 597d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest" 598207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 599d8588912SDave May { 600d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 601d8588912SDave May Vec *L,*R; 602d8588912SDave May MPI_Comm comm; 603d8588912SDave May PetscInt i,j; 604d8588912SDave May PetscErrorCode ierr; 605d8588912SDave May 606d8588912SDave May PetscFunctionBegin; 607d8588912SDave May comm = ((PetscObject)A)->comm; 608d8588912SDave May if (right) { 609d8588912SDave May /* allocate R */ 610d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr); 611d8588912SDave May /* Create the right vectors */ 612d8588912SDave May for (j=0; j<bA->nc; j++) { 613d8588912SDave May for (i=0; i<bA->nr; i++) { 614d8588912SDave May if (bA->m[i][j]) { 615d8588912SDave May ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 616d8588912SDave May break; 617d8588912SDave May } 618d8588912SDave May } 619d8588912SDave May if (i==bA->nr) { 620d8588912SDave May /* have an empty column */ 621d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 622d8588912SDave May } 623d8588912SDave May } 624f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 625d8588912SDave May /* hand back control to the nest vector */ 626d8588912SDave May for (j=0; j<bA->nc; j++) { 6276bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 628d8588912SDave May } 629d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 630d8588912SDave May } 631d8588912SDave May 632d8588912SDave May if (left) { 633d8588912SDave May /* allocate L */ 634d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr); 635d8588912SDave May /* Create the left vectors */ 636d8588912SDave May for (i=0; i<bA->nr; i++) { 637d8588912SDave May for (j=0; j<bA->nc; j++) { 638d8588912SDave May if (bA->m[i][j]) { 639d8588912SDave May ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 640d8588912SDave May break; 641d8588912SDave May } 642d8588912SDave May } 643d8588912SDave May if (j==bA->nc) { 644d8588912SDave May /* have an empty row */ 645d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 646d8588912SDave May } 647d8588912SDave May } 648d8588912SDave May 649f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 650d8588912SDave May for (i=0; i<bA->nr; i++) { 6516bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 652d8588912SDave May } 653d8588912SDave May 654d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 655d8588912SDave May } 656d8588912SDave May PetscFunctionReturn(0); 657d8588912SDave May } 658d8588912SDave May 659d8588912SDave May #undef __FUNCT__ 660d8588912SDave May #define __FUNCT__ "MatView_Nest" 661207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 662d8588912SDave May { 663d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 664d8588912SDave May PetscBool isascii; 665d8588912SDave May PetscInt i,j; 666d8588912SDave May PetscErrorCode ierr; 667d8588912SDave May 668d8588912SDave May PetscFunctionBegin; 669d8588912SDave May ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 670d8588912SDave May if (isascii) { 671d8588912SDave May 672d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n" ); 673d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push0 */ 674d8588912SDave May PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 675d8588912SDave May 676d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" ); 677d8588912SDave May for (i=0; i<bA->nr; i++) { 678d8588912SDave May for (j=0; j<bA->nc; j++) { 679d8588912SDave May const MatType type; 680d8588912SDave May const char *name; 681d8588912SDave May PetscInt NR,NC; 682d8588912SDave May PetscBool isNest = PETSC_FALSE; 683d8588912SDave May 684d8588912SDave May if (!bA->m[i][j]) { 685d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j); 686d8588912SDave May continue; 687d8588912SDave May } 688d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 689d8588912SDave May ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr); 690d8588912SDave May name = ((PetscObject)bA->m[i][j])->prefix; 691d8588912SDave May ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 692d8588912SDave May 693d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC); 694d8588912SDave May 695d8588912SDave May if (isNest) { 696d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push1 */ 697d8588912SDave May ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr); 698d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop1 */ 699d8588912SDave May } 700d8588912SDave May } 701d8588912SDave May } 702d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 703d8588912SDave May } 704d8588912SDave May PetscFunctionReturn(0); 705d8588912SDave May } 706d8588912SDave May 707d8588912SDave May #undef __FUNCT__ 708d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 709207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 710d8588912SDave May { 711d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 712d8588912SDave May PetscInt i,j; 713d8588912SDave May PetscErrorCode ierr; 714d8588912SDave May 715d8588912SDave May PetscFunctionBegin; 716d8588912SDave May for (i=0; i<bA->nr; i++) { 717d8588912SDave May for (j=0; j<bA->nc; j++) { 718d8588912SDave May if (!bA->m[i][j]) continue; 719d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 720d8588912SDave May } 721d8588912SDave May } 722d8588912SDave May PetscFunctionReturn(0); 723d8588912SDave May } 724d8588912SDave May 725d8588912SDave May #undef __FUNCT__ 726d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 727207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 728d8588912SDave May { 729d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 730841e96a3SJed Brown Mat *b; 731841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 732d8588912SDave May PetscErrorCode ierr; 733d8588912SDave May 734d8588912SDave May PetscFunctionBegin; 735841e96a3SJed Brown ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 736841e96a3SJed Brown for (i=0; i<nr; i++) { 737841e96a3SJed Brown for (j=0; j<nc; j++) { 738841e96a3SJed Brown if (bA->m[i][j]) { 739841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 740841e96a3SJed Brown } else { 741841e96a3SJed Brown b[i*nc+j] = PETSC_NULL; 742d8588912SDave May } 743d8588912SDave May } 744d8588912SDave May } 745f349c1fdSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 746841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 747841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 7486bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 749d8588912SDave May } 750d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 751d8588912SDave May 752841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 753841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 754d8588912SDave May PetscFunctionReturn(0); 755d8588912SDave May } 756d8588912SDave May 757d8588912SDave May /* nest api */ 758d8588912SDave May EXTERN_C_BEGIN 759d8588912SDave May #undef __FUNCT__ 760d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 761d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 762d8588912SDave May { 763d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 764d8588912SDave May PetscFunctionBegin; 765d8588912SDave May if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 766d8588912SDave May if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 767d8588912SDave May *mat = bA->m[idxm][jdxm]; 768d8588912SDave May PetscFunctionReturn(0); 769d8588912SDave May } 770d8588912SDave May EXTERN_C_END 771d8588912SDave May 772d8588912SDave May #undef __FUNCT__ 773d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 774d8588912SDave May /*@C 775d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 776d8588912SDave May 777d8588912SDave May Not collective 778d8588912SDave May 779d8588912SDave May Input Parameters: 780629881c0SJed Brown + A - nest matrix 781d8588912SDave May . idxm - index of the matrix within the nest matrix 782629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 783d8588912SDave May 784d8588912SDave May Output Parameter: 785d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 786d8588912SDave May 787d8588912SDave May Level: developer 788d8588912SDave May 789d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 790d8588912SDave May @*/ 7917087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 792d8588912SDave May { 793699a902aSJed Brown PetscErrorCode ierr; 794d8588912SDave May 795d8588912SDave May PetscFunctionBegin; 796699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 797d8588912SDave May PetscFunctionReturn(0); 798d8588912SDave May } 799d8588912SDave May 800d8588912SDave May EXTERN_C_BEGIN 801d8588912SDave May #undef __FUNCT__ 802*0782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 803*0782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 804*0782ca92SJed Brown { 805*0782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 806*0782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 807*0782ca92SJed Brown PetscErrorCode ierr; 808*0782ca92SJed Brown 809*0782ca92SJed Brown PetscFunctionBegin; 810*0782ca92SJed Brown if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 811*0782ca92SJed Brown if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 812*0782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 813*0782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 814*0782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 815*0782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 816*0782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 817*0782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 818*0782ca92SJed Brown if (M != Mi || N != Ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni); 819*0782ca92SJed Brown if (m != mi || n != ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni); 820*0782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 821*0782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 822*0782ca92SJed Brown bA->m[idxm][jdxm] = mat; 823*0782ca92SJed Brown PetscFunctionReturn(0); 824*0782ca92SJed Brown } 825*0782ca92SJed Brown EXTERN_C_END 826*0782ca92SJed Brown 827*0782ca92SJed Brown #undef __FUNCT__ 828*0782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 829*0782ca92SJed Brown /*@C 830*0782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 831*0782ca92SJed Brown 832*0782ca92SJed Brown Logically collective on the submatrix communicator 833*0782ca92SJed Brown 834*0782ca92SJed Brown Input Parameters: 835*0782ca92SJed Brown + A - nest matrix 836*0782ca92SJed Brown . idxm - index of the matrix within the nest matrix 837*0782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 838*0782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 839*0782ca92SJed Brown 840*0782ca92SJed Brown Notes: 841*0782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 842*0782ca92SJed Brown 843*0782ca92SJed Brown This increments the reference count of the submatrix. 844*0782ca92SJed Brown 845*0782ca92SJed Brown Level: developer 846*0782ca92SJed Brown 847*0782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 848*0782ca92SJed Brown @*/ 849*0782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 850*0782ca92SJed Brown { 851*0782ca92SJed Brown PetscErrorCode ierr; 852*0782ca92SJed Brown 853*0782ca92SJed Brown PetscFunctionBegin; 854*0782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 855*0782ca92SJed Brown PetscFunctionReturn(0); 856*0782ca92SJed Brown } 857*0782ca92SJed Brown 858*0782ca92SJed Brown EXTERN_C_BEGIN 859*0782ca92SJed Brown #undef __FUNCT__ 860d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 861d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 862d8588912SDave May { 863d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 864d8588912SDave May PetscFunctionBegin; 865d8588912SDave May if (M) { *M = bA->nr; } 866d8588912SDave May if (N) { *N = bA->nc; } 867d8588912SDave May if (mat) { *mat = bA->m; } 868d8588912SDave May PetscFunctionReturn(0); 869d8588912SDave May } 870d8588912SDave May EXTERN_C_END 871d8588912SDave May 872d8588912SDave May #undef __FUNCT__ 873d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 874d8588912SDave May /*@C 875d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 876d8588912SDave May 877d8588912SDave May Not collective 878d8588912SDave May 879d8588912SDave May Input Parameters: 880629881c0SJed Brown . A - nest matrix 881d8588912SDave May 882d8588912SDave May Output Parameter: 883629881c0SJed Brown + M - number of rows in the nest matrix 884d8588912SDave May . N - number of cols in the nest matrix 885629881c0SJed Brown - mat - 2d array of matrices 886d8588912SDave May 887d8588912SDave May Notes: 888d8588912SDave May 889d8588912SDave May The user should not free the array mat. 890d8588912SDave May 891d8588912SDave May Level: developer 892d8588912SDave May 893d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 894d8588912SDave May @*/ 8957087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 896d8588912SDave May { 897699a902aSJed Brown PetscErrorCode ierr; 898d8588912SDave May 899d8588912SDave May PetscFunctionBegin; 900699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 901d8588912SDave May PetscFunctionReturn(0); 902d8588912SDave May } 903d8588912SDave May 904d8588912SDave May EXTERN_C_BEGIN 905d8588912SDave May #undef __FUNCT__ 906d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 9077087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 908d8588912SDave May { 909d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 910d8588912SDave May 911d8588912SDave May PetscFunctionBegin; 912d8588912SDave May if (M) { *M = bA->nr; } 913d8588912SDave May if (N) { *N = bA->nc; } 914d8588912SDave May PetscFunctionReturn(0); 915d8588912SDave May } 916d8588912SDave May EXTERN_C_END 917d8588912SDave May 918d8588912SDave May #undef __FUNCT__ 919d8588912SDave May #define __FUNCT__ "MatNestGetSize" 920d8588912SDave May /*@C 921d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 922d8588912SDave May 923d8588912SDave May Not collective 924d8588912SDave May 925d8588912SDave May Input Parameters: 926d8588912SDave May . A - nest matrix 927d8588912SDave May 928d8588912SDave May Output Parameter: 929629881c0SJed Brown + M - number of rows in the nested mat 930629881c0SJed Brown - N - number of cols in the nested mat 931d8588912SDave May 932d8588912SDave May Notes: 933d8588912SDave May 934d8588912SDave May Level: developer 935d8588912SDave May 936d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 937d8588912SDave May @*/ 9387087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 939d8588912SDave May { 940699a902aSJed Brown PetscErrorCode ierr; 941d8588912SDave May 942d8588912SDave May PetscFunctionBegin; 943699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 944d8588912SDave May PetscFunctionReturn(0); 945d8588912SDave May } 946d8588912SDave May 947207556f9SJed Brown EXTERN_C_BEGIN 948207556f9SJed Brown #undef __FUNCT__ 949207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 9507087cfbeSBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,const VecType vtype) 951207556f9SJed Brown { 952207556f9SJed Brown PetscErrorCode ierr; 953207556f9SJed Brown PetscBool flg; 954207556f9SJed Brown 955207556f9SJed Brown PetscFunctionBegin; 956207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 957207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 958207556f9SJed Brown A->ops->getvecs = flg ? MatGetVecs_Nest : 0; 9597874fa86SDave May A->ops->getdiagonal = flg ? MatGetDiagonal_Nest : 0; 9607874fa86SDave May A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0; 9617874fa86SDave May 962207556f9SJed Brown PetscFunctionReturn(0); 963207556f9SJed Brown } 964207556f9SJed Brown EXTERN_C_END 965207556f9SJed Brown 966207556f9SJed Brown #undef __FUNCT__ 967207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 968207556f9SJed Brown /*@C 969207556f9SJed Brown MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 970207556f9SJed Brown 971207556f9SJed Brown Not collective 972207556f9SJed Brown 973207556f9SJed Brown Input Parameters: 974207556f9SJed Brown + A - nest matrix 975207556f9SJed Brown - vtype - type to use for creating vectors 976207556f9SJed Brown 977207556f9SJed Brown Notes: 978207556f9SJed Brown 979207556f9SJed Brown Level: developer 980207556f9SJed Brown 981207556f9SJed Brown .seealso: MatGetVecs() 982207556f9SJed Brown @*/ 9837087cfbeSBarry Smith PetscErrorCode MatNestSetVecType(Mat A,const VecType vtype) 984207556f9SJed Brown { 985207556f9SJed Brown PetscErrorCode ierr; 986207556f9SJed Brown 987207556f9SJed Brown PetscFunctionBegin; 988207556f9SJed Brown ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr); 989207556f9SJed Brown PetscFunctionReturn(0); 990207556f9SJed Brown } 991207556f9SJed Brown 992c8883902SJed Brown EXTERN_C_BEGIN 993d8588912SDave May #undef __FUNCT__ 994c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 995c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 996d8588912SDave May { 997c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 998c8883902SJed Brown PetscInt i,j,m,n,M,N; 999d8588912SDave May PetscErrorCode ierr; 1000d8588912SDave May 1001d8588912SDave May PetscFunctionBegin; 1002c8883902SJed Brown s->nr = nr; 1003c8883902SJed Brown s->nc = nc; 1004d8588912SDave May 1005c8883902SJed Brown /* Create space for submatrices */ 1006c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 1007c8883902SJed Brown for (i=0; i<nr; i++) { 1008c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 1009d8588912SDave May } 1010c8883902SJed Brown for (i=0; i<nr; i++) { 1011c8883902SJed Brown for (j=0; j<nc; j++) { 1012c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1013c8883902SJed Brown if (a[i*nc+j]) { 1014c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1015d8588912SDave May } 1016d8588912SDave May } 1017d8588912SDave May } 1018d8588912SDave May 10198188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1020d8588912SDave May 1021c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 1022c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 1023c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1024c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1025d8588912SDave May 10268188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1027d8588912SDave May 1028c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1029c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1030c8883902SJed Brown ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 1031c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1032c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1033c8883902SJed Brown ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 1034c8883902SJed Brown 1035c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1036c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1037c8883902SJed Brown 1038c8883902SJed Brown ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 1039c8883902SJed Brown ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr); 1040c8883902SJed Brown ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr); 1041d8588912SDave May PetscFunctionReturn(0); 1042d8588912SDave May } 1043c8883902SJed Brown EXTERN_C_END 1044d8588912SDave May 1045c8883902SJed Brown #undef __FUNCT__ 1046c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 1047c8883902SJed Brown /*@ 1048c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1049c8883902SJed Brown 1050c8883902SJed Brown Collective on Mat 1051c8883902SJed Brown 1052c8883902SJed Brown Input Parameter: 1053c8883902SJed Brown + N - nested matrix 1054c8883902SJed Brown . nr - number of nested row blocks 1055c8883902SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1056c8883902SJed Brown . nc - number of nested column blocks 1057c8883902SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1058c8883902SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1059c8883902SJed Brown 1060c8883902SJed Brown Level: advanced 1061c8883902SJed Brown 1062c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 1063c8883902SJed Brown @*/ 1064c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1065c8883902SJed Brown { 1066c8883902SJed Brown PetscErrorCode ierr; 1067c8883902SJed Brown PetscInt i; 1068c8883902SJed Brown 1069c8883902SJed Brown PetscFunctionBegin; 1070c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1071c8883902SJed Brown if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1072c8883902SJed Brown if (nr && is_row) { 1073c8883902SJed Brown PetscValidPointer(is_row,3); 1074c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1075c8883902SJed Brown } 1076c8883902SJed Brown if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 1077c8883902SJed Brown if (nc && is_row) { 1078c8883902SJed Brown PetscValidPointer(is_col,5); 1079c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1080c8883902SJed Brown } 1081c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 1082c8883902SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr); 1083c8883902SJed Brown PetscFunctionReturn(0); 1084c8883902SJed Brown } 1085d8588912SDave May 1086d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1087d8588912SDave May /* 1088d8588912SDave May nprocessors = NP 1089d8588912SDave May Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) ) 1090d8588912SDave May proc 0: => (g_0,h_0,) 1091d8588912SDave May proc 1: => (g_1,h_1,) 1092d8588912SDave May ... 1093d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1094d8588912SDave May 1095d8588912SDave May proc 0: proc 1: proc nprocs-1: 1096d8588912SDave May is[0] = ( 0,1,2,...,nlocal(g_0)-1 ) ( 0,1,...,nlocal(g_1)-1 ) ( 0,1,...,nlocal(g_NP-1) ) 1097d8588912SDave May 1098d8588912SDave May proc 0: 1099d8588912SDave May is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 ) 1100d8588912SDave May proc 1: 1101d8588912SDave May is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 ) 1102d8588912SDave May 1103d8588912SDave May proc NP-1: 1104d8588912SDave May is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 ) 1105d8588912SDave May */ 1106d8588912SDave May #undef __FUNCT__ 1107d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1108841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1109d8588912SDave May { 1110e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 11118188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1112d8588912SDave May PetscErrorCode ierr; 11132ae74bdbSJed Brown Mat sub; 1114d8588912SDave May 1115d8588912SDave May PetscFunctionBegin; 11168188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 11178188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1118d8588912SDave May if (is_row) { /* valid IS is passed in */ 1119d8588912SDave May /* refs on is[] are incremeneted */ 1120e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1121d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1122e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1123d8588912SDave May } 11242ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 11258188e55aSJed Brown nsum = 0; 11268188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 11278188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 11288188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 11298188e55aSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 11308188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 11318188e55aSJed Brown nsum += n; 11328188e55aSJed Brown } 11338188e55aSJed Brown offset = 0; 11348188e55aSJed Brown ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1135e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1136f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 11372ae74bdbSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 11382ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1139e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1140e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 11412ae74bdbSJed Brown offset += n; 1142d8588912SDave May } 1143d8588912SDave May } 1144d8588912SDave May 1145d8588912SDave May if (is_col) { /* valid IS is passed in */ 1146d8588912SDave May /* refs on is[] are incremeneted */ 1147e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1148d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1149e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1150d8588912SDave May } 11512ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 11522ae74bdbSJed Brown offset = A->cmap->rstart; 11538188e55aSJed Brown nsum = 0; 11548188e55aSJed Brown for (j=0; j<vs->nc; j++) { 11558188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 11568188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 11578188e55aSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 11588188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 11598188e55aSJed Brown nsum += n; 11608188e55aSJed Brown } 11618188e55aSJed Brown offset = 0; 11628188e55aSJed Brown ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1163e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1164f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 11652ae74bdbSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 11662ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1167e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1168e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 11692ae74bdbSJed Brown offset += n; 1170d8588912SDave May } 1171d8588912SDave May } 1172e2d7f03fSJed Brown 1173e2d7f03fSJed Brown /* Set up the local ISs */ 1174e2d7f03fSJed Brown ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1175e2d7f03fSJed Brown ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1176e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1177e2d7f03fSJed Brown IS isloc; 11788188e55aSJed Brown ISLocalToGlobalMapping rmap = PETSC_NULL; 1179e2d7f03fSJed Brown PetscInt nlocal,bs; 1180e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 11818188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);} 1182207556f9SJed Brown if (rmap) { 1183e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1184e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1185e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1186e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1187207556f9SJed Brown } else { 1188207556f9SJed Brown nlocal = 0; 1189207556f9SJed Brown isloc = PETSC_NULL; 1190207556f9SJed Brown } 1191e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1192e2d7f03fSJed Brown offset += nlocal; 1193e2d7f03fSJed Brown } 11948188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1195e2d7f03fSJed Brown IS isloc; 11968188e55aSJed Brown ISLocalToGlobalMapping cmap = PETSC_NULL; 1197e2d7f03fSJed Brown PetscInt nlocal,bs; 1198e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 11998188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);} 1200207556f9SJed Brown if (cmap) { 1201e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1202e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1203e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1204e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1205207556f9SJed Brown } else { 1206207556f9SJed Brown nlocal = 0; 1207207556f9SJed Brown isloc = PETSC_NULL; 1208207556f9SJed Brown } 1209e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1210e2d7f03fSJed Brown offset += nlocal; 1211e2d7f03fSJed Brown } 1212d8588912SDave May PetscFunctionReturn(0); 1213d8588912SDave May } 1214d8588912SDave May 1215d8588912SDave May #undef __FUNCT__ 1216d8588912SDave May #define __FUNCT__ "MatCreateNest" 1217659c6bb0SJed Brown /*@ 1218659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1219659c6bb0SJed Brown 1220659c6bb0SJed Brown Collective on Mat 1221659c6bb0SJed Brown 1222659c6bb0SJed Brown Input Parameter: 1223659c6bb0SJed Brown + comm - Communicator for the new Mat 1224659c6bb0SJed Brown . nr - number of nested row blocks 1225659c6bb0SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1226659c6bb0SJed Brown . nc - number of nested column blocks 1227659c6bb0SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1228659c6bb0SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1229659c6bb0SJed Brown 1230659c6bb0SJed Brown Output Parameter: 1231659c6bb0SJed Brown . B - new matrix 1232659c6bb0SJed Brown 1233659c6bb0SJed Brown Level: advanced 1234659c6bb0SJed Brown 1235659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST 1236659c6bb0SJed Brown @*/ 12377087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1238d8588912SDave May { 1239d8588912SDave May Mat A; 1240d8588912SDave May PetscErrorCode ierr; 1241d8588912SDave May 1242d8588912SDave May PetscFunctionBegin; 1243c8883902SJed Brown *B = 0; 1244d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1245c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 1246c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1247d8588912SDave May *B = A; 1248d8588912SDave May PetscFunctionReturn(0); 1249d8588912SDave May } 1250659c6bb0SJed Brown 1251659c6bb0SJed Brown /*MC 1252659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1253659c6bb0SJed Brown 1254659c6bb0SJed Brown Level: intermediate 1255659c6bb0SJed Brown 1256659c6bb0SJed Brown Notes: 1257659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1258659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1259659c6bb0SJed Brown It is usually used with DMComposite and DMGetMatrix() 1260659c6bb0SJed Brown 1261659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1262659c6bb0SJed Brown M*/ 1263c8883902SJed Brown EXTERN_C_BEGIN 1264c8883902SJed Brown #undef __FUNCT__ 1265c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 1266c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A) 1267c8883902SJed Brown { 1268c8883902SJed Brown Mat_Nest *s; 1269c8883902SJed Brown PetscErrorCode ierr; 1270c8883902SJed Brown 1271c8883902SJed Brown PetscFunctionBegin; 1272c8883902SJed Brown ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr); 1273c8883902SJed Brown A->data = (void*)s; 1274c8883902SJed Brown s->nr = s->nc = -1; 1275c8883902SJed Brown s->m = PETSC_NULL; 1276c8883902SJed Brown 1277c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 1278c8883902SJed Brown A->ops->mult = MatMult_Nest; 1279c8883902SJed Brown A->ops->multadd = 0; 1280c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 1281c8883902SJed Brown A->ops->multtransposeadd = 0; 1282c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1283c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1284c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1285c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1286c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1287c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1288c8883902SJed Brown A->ops->view = MatView_Nest; 1289c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1290c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1291c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 12927874fa86SDave May A->ops->getdiagonal = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */ 12937874fa86SDave May A->ops->diagonalscale = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */ 1294c8883902SJed Brown 1295c8883902SJed Brown A->spptr = 0; 1296c8883902SJed Brown A->same_nonzero = PETSC_FALSE; 1297c8883902SJed Brown A->assembled = PETSC_FALSE; 1298c8883902SJed Brown 1299c8883902SJed Brown /* expose Nest api's */ 1300c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "MatNestGetSubMat_Nest", MatNestGetSubMat_Nest);CHKERRQ(ierr); 1301*0782ca92SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "MatNestSetSubMat_Nest", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1302c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "MatNestGetSubMats_Nest", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1303c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "MatNestGetSize_Nest", MatNestGetSize_Nest);CHKERRQ(ierr); 1304c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "MatNestSetVecType_Nest", MatNestSetVecType_Nest);CHKERRQ(ierr); 1305c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "MatNestSetSubMats_Nest", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1306c8883902SJed Brown 1307c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1308c8883902SJed Brown PetscFunctionReturn(0); 1309c8883902SJed Brown } 131086e80854SHong Zhang EXTERN_C_END 1311