1d8588912SDave May #define PETSCMAT_DLL 2d8588912SDave May 3ec9811eeSJed Brown #include "matnestimpl.h" /*I "petscmat.h" I*/ 4d8588912SDave May 5d8588912SDave May /* private functions */ 6d8588912SDave May #undef __FUNCT__ 7d8588912SDave May #define __FUNCT__ "MatNestGetSize_Recursive" 8d8588912SDave May static PetscErrorCode MatNestGetSize_Recursive(Mat A,PetscBool globalsize,PetscInt *M,PetscInt *N) 9d8588912SDave May { 10d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 11*b9dfa011SHong Zhang PetscInt sizeM,sizeN,i,j,index=-1; 12d8588912SDave May PetscErrorCode ierr; 13d8588912SDave May 14d8588912SDave May PetscFunctionBegin; 15d8588912SDave May /* rows */ 16d8588912SDave May /* now descend recursively */ 17d8588912SDave May for (i=0; i<bA->nr; i++) { 18d8588912SDave May for (j=0; j<bA->nc; j++) { 19d8588912SDave May if (bA->m[i][j]) { index = j; break; } 20d8588912SDave May } 21d8588912SDave May if (bA->m[i][index]) { 22d8588912SDave May if (globalsize) { ierr = MatGetSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); } 23d8588912SDave May else { ierr = MatGetLocalSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); } 24d8588912SDave May *M = *M + sizeM; 25d8588912SDave May } 26d8588912SDave May } 27d8588912SDave May 28d8588912SDave May /* cols */ 29d8588912SDave May /* now descend recursively */ 30d8588912SDave May for (j=0; j<bA->nc; j++) { 31d8588912SDave May for (i=0; i<bA->nr; i++) { 32d8588912SDave May if (bA->m[i][j]) { index = i; break; } 33d8588912SDave May } 34d8588912SDave May if (bA->m[index][j]) { 35d8588912SDave May if (globalsize) { ierr = MatGetSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); } 36d8588912SDave May else { ierr = MatGetLocalSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); } 37d8588912SDave May *N = *N + sizeN; 38d8588912SDave May } 39d8588912SDave May } 40d8588912SDave May PetscFunctionReturn(0); 41d8588912SDave May } 42d8588912SDave May 43d8588912SDave May #undef __FUNCT__ 44d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2" 45207556f9SJed Brown PETSC_UNUSED 46d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y) 47d8588912SDave May { 48d8588912SDave May PetscBool isANest, isxNest, isyNest; 49d8588912SDave May PetscErrorCode ierr; 50d8588912SDave May 51d8588912SDave May PetscFunctionBegin; 52d8588912SDave May isANest = isxNest = isyNest = PETSC_FALSE; 53d8588912SDave May ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr); 54d8588912SDave May ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr); 55d8588912SDave May ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr); 56d8588912SDave May 57d8588912SDave May if (isANest && isxNest && isyNest) { 58d8588912SDave May PetscFunctionReturn(0); 59d8588912SDave May } else { 60d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations"); 61d8588912SDave May PetscFunctionReturn(0); 62d8588912SDave May } 63d8588912SDave May PetscFunctionReturn(0); 64d8588912SDave May } 65d8588912SDave May 66d8588912SDave May /* 67d8588912SDave May This function is called every time we insert a sub matrix the Nest. 68d8588912SDave May It ensures that every Nest along row r, and col c has the same dimensions 69d8588912SDave May as the submat being inserted. 70d8588912SDave May */ 71d8588912SDave May #undef __FUNCT__ 72d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckConsistency" 73063a28d4SJed Brown PETSC_UNUSED 74d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c) 75d8588912SDave May { 76d8588912SDave May Mat_Nest *b = (Mat_Nest*)A->data; 77d8588912SDave May PetscInt i,j; 78d8588912SDave May PetscInt nr,nc; 79d8588912SDave May PetscInt sM,sN,M,N; 80d8588912SDave May Mat mat; 81d8588912SDave May PetscErrorCode ierr; 82d8588912SDave May 83d8588912SDave May PetscFunctionBegin; 84d8588912SDave May nr = b->nr; 85d8588912SDave May nc = b->nc; 86d8588912SDave May ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr); 87d8588912SDave May for (i=0; i<nr; i++) { 88d8588912SDave May mat = b->m[i][c]; 89d8588912SDave May if (mat) { 90d8588912SDave May ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 91d8588912SDave May /* Check columns match */ 92d8588912SDave May if (sN != N) { 93d8588912SDave May SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N ); 94d8588912SDave May } 95d8588912SDave May } 96d8588912SDave May } 97d8588912SDave May 98d8588912SDave May for (j=0; j<nc; j++) { 99d8588912SDave May mat = b->m[r][j]; 100d8588912SDave May if (mat) { 101d8588912SDave May ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 102d8588912SDave May /* Check rows match */ 103d8588912SDave May if (sM != M) { 104d8588912SDave May SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M ); 105d8588912SDave May } 106d8588912SDave May } 107d8588912SDave May } 108d8588912SDave May PetscFunctionReturn(0); 109d8588912SDave May } 110d8588912SDave May 111d8588912SDave May /* 112d8588912SDave May Checks if the row/col's contain a complete set of IS's. 113d8588912SDave May Once they are filled, the offsets are computed. This saves having to update 114d8588912SDave May the index set entries each time we insert something new. 115d8588912SDave May */ 116d8588912SDave May #undef __FUNCT__ 117d8588912SDave May #define __FUNCT__ "PETSc_MatNest_UpdateStructure" 118063a28d4SJed Brown PETSC_UNUSED 119d8588912SDave May static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx) 120d8588912SDave May { 121d8588912SDave May Mat_Nest *b = (Mat_Nest*)A->data; 122d8588912SDave May PetscInt m,n,i,j; 123d8588912SDave May PetscBool fullrow,fullcol; 124d8588912SDave May PetscErrorCode ierr; 125d8588912SDave May 126d8588912SDave May PetscFunctionBegin; 127d8588912SDave May ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr); 128d8588912SDave May b->row_len[ridx] = m; 129d8588912SDave May b->col_len[cidx] = n; 130d8588912SDave May 131d8588912SDave May fullrow = PETSC_TRUE; 132d8588912SDave May for (i=0; i<b->nr; i++) { 133d8588912SDave May if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; } 134d8588912SDave May } 135f349c1fdSJed Brown if ( (fullrow) && (!b->isglobal.row[0]) ){ 136d8588912SDave May PetscInt cnt = 0; 137d8588912SDave May for (i=0; i<b->nr; i++) { 138f349c1fdSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr); 139d8588912SDave May cnt = cnt + b->row_len[i]; 140d8588912SDave May } 141d8588912SDave May } 142d8588912SDave May 143d8588912SDave May fullcol = PETSC_TRUE; 144d8588912SDave May for (j=0; j<b->nc; j++) { 145d8588912SDave May if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; } 146d8588912SDave May } 147f349c1fdSJed Brown if( (fullcol) && (!b->isglobal.col[0]) ){ 148d8588912SDave May PetscInt cnt = 0; 149d8588912SDave May for (j=0; j<b->nc; j++) { 150f349c1fdSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr); 151d8588912SDave May cnt = cnt + b->col_len[j]; 152d8588912SDave May } 153d8588912SDave May } 154d8588912SDave May PetscFunctionReturn(0); 155d8588912SDave May } 156d8588912SDave May 157d8588912SDave May /* operations */ 158d8588912SDave May #undef __FUNCT__ 159d8588912SDave May #define __FUNCT__ "MatMult_Nest" 160207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 161d8588912SDave May { 162d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 163207556f9SJed Brown Vec *bx = bA->right,*by = bA->left; 164207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 165d8588912SDave May PetscErrorCode ierr; 166d8588912SDave May 167d8588912SDave May PetscFunctionBegin; 168207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 169207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 170207556f9SJed Brown for (i=0; i<nr; i++) { 171d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 172207556f9SJed Brown for (j=0; j<nc; j++) { 173207556f9SJed Brown if (!bA->m[i][j]) continue; 174d8588912SDave May /* y[i] <- y[i] + A[i][j] * x[j] */ 175d8588912SDave May ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 176d8588912SDave May } 177d8588912SDave May } 178207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 179207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 180d8588912SDave May PetscFunctionReturn(0); 181d8588912SDave May } 182d8588912SDave May 183d8588912SDave May #undef __FUNCT__ 184d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest" 185207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 186d8588912SDave May { 187d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 188207556f9SJed Brown Vec *bx = bA->left,*by = bA->right; 189207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 190d8588912SDave May PetscErrorCode ierr; 191d8588912SDave May 192d8588912SDave May PetscFunctionBegin; 193d8588912SDave May if (A->symmetric) { 194d8588912SDave May ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr); 195d8588912SDave May PetscFunctionReturn(0); 196d8588912SDave May } 197207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 198207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 199207556f9SJed Brown for (i=0; i<nr; i++) { 200d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 201207556f9SJed Brown for (j=0; j<nc; j++) { 202207556f9SJed Brown if (!bA->m[j][i]) continue; 203d8588912SDave May /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */ 204d8588912SDave May ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr); 205d8588912SDave May } 206d8588912SDave May } 207207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 208207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 209d8588912SDave May PetscFunctionReturn(0); 210d8588912SDave May } 211d8588912SDave May 212d8588912SDave May #undef __FUNCT__ 213e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList" 214e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 215e2d7f03fSJed Brown { 216e2d7f03fSJed Brown PetscErrorCode ierr; 217e2d7f03fSJed Brown IS *lst = *list; 218e2d7f03fSJed Brown PetscInt i; 219e2d7f03fSJed Brown 220e2d7f03fSJed Brown PetscFunctionBegin; 221e2d7f03fSJed Brown if (!lst) PetscFunctionReturn(0); 222e2d7f03fSJed Brown for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(lst[i]);CHKERRQ(ierr);} 223e2d7f03fSJed Brown ierr = PetscFree(lst);CHKERRQ(ierr); 224e2d7f03fSJed Brown *list = PETSC_NULL; 225e2d7f03fSJed Brown PetscFunctionReturn(0); 226e2d7f03fSJed Brown } 227e2d7f03fSJed Brown 228e2d7f03fSJed Brown #undef __FUNCT__ 229d8588912SDave May #define __FUNCT__ "MatDestroy_Nest" 230207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A) 231d8588912SDave May { 232d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 233d8588912SDave May PetscInt i,j; 234d8588912SDave May PetscErrorCode ierr; 235d8588912SDave May 236d8588912SDave May PetscFunctionBegin; 237d8588912SDave May /* release the matrices and the place holders */ 238e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 239e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 240e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 241e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 242d8588912SDave May 243d8588912SDave May ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 244d8588912SDave May ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 245d8588912SDave May 246207556f9SJed Brown ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 247207556f9SJed Brown 248d8588912SDave May /* release the matrices and the place holders */ 249d8588912SDave May if (vs->m) { 250d8588912SDave May for (i=0; i<vs->nr; i++) { 251d8588912SDave May for (j=0; j<vs->nc; j++) { 252d8588912SDave May 253d8588912SDave May if (vs->m[i][j]) { 254d8588912SDave May ierr = MatDestroy(vs->m[i][j]);CHKERRQ(ierr); 255d8588912SDave May vs->m[i][j] = PETSC_NULL; 256d8588912SDave May } 257d8588912SDave May } 258d8588912SDave May ierr = PetscFree( vs->m[i] );CHKERRQ(ierr); 259d8588912SDave May vs->m[i] = PETSC_NULL; 260d8588912SDave May } 261d8588912SDave May ierr = PetscFree(vs->m);CHKERRQ(ierr); 262d8588912SDave May vs->m = PETSC_NULL; 263d8588912SDave May } 264d8588912SDave May ierr = PetscFree(vs);CHKERRQ(ierr); 265d8588912SDave May 266d8588912SDave May PetscFunctionReturn(0); 267d8588912SDave May } 268d8588912SDave May 269d8588912SDave May #undef __FUNCT__ 270d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest" 271207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 272d8588912SDave May { 273d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 274d8588912SDave May PetscInt i,j; 275d8588912SDave May PetscErrorCode ierr; 276d8588912SDave May 277d8588912SDave May PetscFunctionBegin; 278d8588912SDave May for (i=0; i<vs->nr; i++) { 279d8588912SDave May for (j=0; j<vs->nc; j++) { 280d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); } 281d8588912SDave May } 282d8588912SDave May } 283d8588912SDave May PetscFunctionReturn(0); 284d8588912SDave May } 285d8588912SDave May 286d8588912SDave May #undef __FUNCT__ 287d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 288207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 289d8588912SDave May { 290d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 291d8588912SDave May PetscInt i,j; 292d8588912SDave May PetscErrorCode ierr; 293d8588912SDave May 294d8588912SDave May PetscFunctionBegin; 295d8588912SDave May for (i=0; i<vs->nr; i++) { 296d8588912SDave May for (j=0; j<vs->nc; j++) { 297d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); } 298d8588912SDave May } 299d8588912SDave May } 300d8588912SDave May PetscFunctionReturn(0); 301d8588912SDave May } 302d8588912SDave May 303d8588912SDave May #undef __FUNCT__ 304f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 305f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 306d8588912SDave May { 307207556f9SJed Brown PetscErrorCode ierr; 308f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 309f349c1fdSJed Brown PetscInt j; 310f349c1fdSJed Brown Mat sub; 311d8588912SDave May 312d8588912SDave May PetscFunctionBegin; 313f349c1fdSJed Brown sub = vs->m[row][row]; /* Prefer to find on the diagonal */ 314f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 315f349c1fdSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",row); 316207556f9SJed Brown ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */ 317f349c1fdSJed Brown *B = sub; 318f349c1fdSJed Brown PetscFunctionReturn(0); 319d8588912SDave May } 320d8588912SDave May 321f349c1fdSJed Brown #undef __FUNCT__ 322f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 323f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 324f349c1fdSJed Brown { 325207556f9SJed Brown PetscErrorCode ierr; 326f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 327f349c1fdSJed Brown PetscInt i; 328f349c1fdSJed Brown Mat sub; 329f349c1fdSJed Brown 330f349c1fdSJed Brown PetscFunctionBegin; 331f349c1fdSJed Brown sub = vs->m[col][col]; /* Prefer to find on the diagonal */ 332f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 333f349c1fdSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",col); 334207556f9SJed Brown ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */ 335f349c1fdSJed Brown *B = sub; 336f349c1fdSJed Brown PetscFunctionReturn(0); 337d8588912SDave May } 338d8588912SDave May 339f349c1fdSJed Brown #undef __FUNCT__ 340f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 341f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 342f349c1fdSJed Brown { 343f349c1fdSJed Brown PetscErrorCode ierr; 344f349c1fdSJed Brown PetscInt i; 345f349c1fdSJed Brown PetscBool flg; 346f349c1fdSJed Brown 347f349c1fdSJed Brown PetscFunctionBegin; 348f349c1fdSJed Brown PetscValidPointer(list,3); 349f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 350f349c1fdSJed Brown PetscValidIntPointer(found,5); 351f349c1fdSJed Brown *found = -1; 352f349c1fdSJed Brown for (i=0; i<n; i++) { 353207556f9SJed Brown if (!list[i]) continue; 354f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 355f349c1fdSJed Brown if (flg) { 356f349c1fdSJed Brown *found = i; 357f349c1fdSJed Brown PetscFunctionReturn(0); 358f349c1fdSJed Brown } 359f349c1fdSJed Brown } 360f349c1fdSJed Brown SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set"); 361f349c1fdSJed Brown PetscFunctionReturn(0); 362f349c1fdSJed Brown } 363f349c1fdSJed Brown 364f349c1fdSJed Brown #undef __FUNCT__ 365f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 366f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 367f349c1fdSJed Brown { 368f349c1fdSJed Brown PetscErrorCode ierr; 369f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 370f349c1fdSJed Brown PetscInt row,col; 371f349c1fdSJed Brown 372f349c1fdSJed Brown PetscFunctionBegin; 373f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 374f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 375f349c1fdSJed Brown *B = vs->m[row][col]; 376f349c1fdSJed Brown PetscFunctionReturn(0); 377f349c1fdSJed Brown } 378f349c1fdSJed Brown 379f349c1fdSJed Brown #undef __FUNCT__ 380f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 381f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 382f349c1fdSJed Brown { 383f349c1fdSJed Brown PetscErrorCode ierr; 384f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 385f349c1fdSJed Brown Mat sub; 386f349c1fdSJed Brown 387f349c1fdSJed Brown PetscFunctionBegin; 388f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 389f349c1fdSJed Brown switch (reuse) { 390f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 391f349c1fdSJed Brown ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); 392f349c1fdSJed Brown *B = sub; 393f349c1fdSJed Brown break; 394f349c1fdSJed Brown case MAT_REUSE_MATRIX: 395f349c1fdSJed Brown if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 396f349c1fdSJed Brown break; 397f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 398f349c1fdSJed Brown break; 399f349c1fdSJed Brown } 400f349c1fdSJed Brown PetscFunctionReturn(0); 401f349c1fdSJed Brown } 402f349c1fdSJed Brown 403f349c1fdSJed Brown #undef __FUNCT__ 404f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 405f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 406f349c1fdSJed Brown { 407f349c1fdSJed Brown PetscErrorCode ierr; 408f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 409f349c1fdSJed Brown Mat sub; 410f349c1fdSJed Brown 411f349c1fdSJed Brown PetscFunctionBegin; 412f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 413f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 414f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 415f349c1fdSJed Brown *B = sub; 416d8588912SDave May PetscFunctionReturn(0); 417d8588912SDave May } 418d8588912SDave May 419d8588912SDave May #undef __FUNCT__ 420d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 421207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 422d8588912SDave May { 423d8588912SDave May PetscErrorCode ierr; 424f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 425f349c1fdSJed Brown Mat sub; 426d8588912SDave May 427d8588912SDave May PetscFunctionBegin; 428f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 429f349c1fdSJed Brown if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 430f349c1fdSJed Brown if (sub) { 431f349c1fdSJed Brown if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 432f349c1fdSJed Brown ierr = MatDestroy(*B);CHKERRQ(ierr); 433d8588912SDave May } 434f349c1fdSJed Brown *B = PETSC_NULL; 435d8588912SDave May PetscFunctionReturn(0); 436d8588912SDave May } 437d8588912SDave May 438d8588912SDave May #undef __FUNCT__ 439d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest" 440207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 441d8588912SDave May { 442d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 443d8588912SDave May Vec *L,*R; 444d8588912SDave May MPI_Comm comm; 445d8588912SDave May PetscInt i,j; 446d8588912SDave May PetscErrorCode ierr; 447d8588912SDave May 448d8588912SDave May PetscFunctionBegin; 449d8588912SDave May comm = ((PetscObject)A)->comm; 450d8588912SDave May if (right) { 451d8588912SDave May /* allocate R */ 452d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr); 453d8588912SDave May /* Create the right vectors */ 454d8588912SDave May for (j=0; j<bA->nc; j++) { 455d8588912SDave May for (i=0; i<bA->nr; i++) { 456d8588912SDave May if (bA->m[i][j]) { 457d8588912SDave May ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 458d8588912SDave May break; 459d8588912SDave May } 460d8588912SDave May } 461d8588912SDave May if (i==bA->nr) { 462d8588912SDave May /* have an empty column */ 463d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 464d8588912SDave May } 465d8588912SDave May } 466f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 467d8588912SDave May /* hand back control to the nest vector */ 468d8588912SDave May for (j=0; j<bA->nc; j++) { 469d8588912SDave May ierr = VecDestroy(R[j]);CHKERRQ(ierr); 470d8588912SDave May } 471d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 472d8588912SDave May } 473d8588912SDave May 474d8588912SDave May if (left) { 475d8588912SDave May /* allocate L */ 476d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr); 477d8588912SDave May /* Create the left vectors */ 478d8588912SDave May for (i=0; i<bA->nr; i++) { 479d8588912SDave May for (j=0; j<bA->nc; j++) { 480d8588912SDave May if (bA->m[i][j]) { 481d8588912SDave May ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 482d8588912SDave May break; 483d8588912SDave May } 484d8588912SDave May } 485d8588912SDave May if (j==bA->nc) { 486d8588912SDave May /* have an empty row */ 487d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 488d8588912SDave May } 489d8588912SDave May } 490d8588912SDave May 491f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 492d8588912SDave May for (i=0; i<bA->nr; i++) { 493d8588912SDave May ierr = VecDestroy(L[i]);CHKERRQ(ierr); 494d8588912SDave May } 495d8588912SDave May 496d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 497d8588912SDave May } 498d8588912SDave May PetscFunctionReturn(0); 499d8588912SDave May } 500d8588912SDave May 501d8588912SDave May #undef __FUNCT__ 502d8588912SDave May #define __FUNCT__ "MatView_Nest" 503207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 504d8588912SDave May { 505d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 506d8588912SDave May PetscBool isascii; 507d8588912SDave May PetscInt i,j; 508d8588912SDave May PetscErrorCode ierr; 509d8588912SDave May 510d8588912SDave May PetscFunctionBegin; 511d8588912SDave May ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 512d8588912SDave May if (isascii) { 513d8588912SDave May 514d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n" ); 515d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push0 */ 516d8588912SDave May PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 517d8588912SDave May 518d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" ); 519d8588912SDave May for (i=0; i<bA->nr; i++) { 520d8588912SDave May for (j=0; j<bA->nc; j++) { 521d8588912SDave May const MatType type; 522d8588912SDave May const char *name; 523d8588912SDave May PetscInt NR,NC; 524d8588912SDave May PetscBool isNest = PETSC_FALSE; 525d8588912SDave May 526d8588912SDave May if (!bA->m[i][j]) { 527d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j); 528d8588912SDave May continue; 529d8588912SDave May } 530d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 531d8588912SDave May ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr); 532d8588912SDave May name = ((PetscObject)bA->m[i][j])->prefix; 533d8588912SDave May ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 534d8588912SDave May 535d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC); 536d8588912SDave May 537d8588912SDave May if (isNest) { 538d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push1 */ 539d8588912SDave May ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr); 540d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop1 */ 541d8588912SDave May } 542d8588912SDave May } 543d8588912SDave May } 544d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 545d8588912SDave May } 546d8588912SDave May PetscFunctionReturn(0); 547d8588912SDave May } 548d8588912SDave May 549d8588912SDave May #undef __FUNCT__ 550d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 551207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 552d8588912SDave May { 553d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 554d8588912SDave May PetscInt i,j; 555d8588912SDave May PetscErrorCode ierr; 556d8588912SDave May 557d8588912SDave May PetscFunctionBegin; 558d8588912SDave May for (i=0; i<bA->nr; i++) { 559d8588912SDave May for (j=0; j<bA->nc; j++) { 560d8588912SDave May if (!bA->m[i][j]) continue; 561d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 562d8588912SDave May } 563d8588912SDave May } 564d8588912SDave May PetscFunctionReturn(0); 565d8588912SDave May } 566d8588912SDave May 567d8588912SDave May #undef __FUNCT__ 568d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 569207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 570d8588912SDave May { 571d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 572841e96a3SJed Brown Mat *b; 573841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 574d8588912SDave May PetscErrorCode ierr; 575d8588912SDave May 576d8588912SDave May PetscFunctionBegin; 577841e96a3SJed Brown ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 578841e96a3SJed Brown for (i=0; i<nr; i++) { 579841e96a3SJed Brown for (j=0; j<nc; j++) { 580841e96a3SJed Brown if (bA->m[i][j]) { 581841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 582841e96a3SJed Brown } else { 583841e96a3SJed Brown b[i*nc+j] = PETSC_NULL; 584d8588912SDave May } 585d8588912SDave May } 586d8588912SDave May } 587f349c1fdSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 588841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 589841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 590841e96a3SJed Brown if (b[i]) {ierr = MatDestroy(b[i]);CHKERRQ(ierr);} 591d8588912SDave May } 592d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 593d8588912SDave May 594841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 595841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 596d8588912SDave May PetscFunctionReturn(0); 597d8588912SDave May } 598d8588912SDave May 599d8588912SDave May /* nest api */ 600d8588912SDave May EXTERN_C_BEGIN 601d8588912SDave May #undef __FUNCT__ 602d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 603d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 604d8588912SDave May { 605d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 606d8588912SDave May PetscFunctionBegin; 607d8588912SDave May if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 608d8588912SDave May if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 609d8588912SDave May *mat = bA->m[idxm][jdxm]; 610d8588912SDave May PetscFunctionReturn(0); 611d8588912SDave May } 612d8588912SDave May EXTERN_C_END 613d8588912SDave May 614d8588912SDave May #undef __FUNCT__ 615d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 616d8588912SDave May /*@C 617d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 618d8588912SDave May 619d8588912SDave May Not collective 620d8588912SDave May 621d8588912SDave May Input Parameters: 622d8588912SDave May . A - nest matrix 623d8588912SDave May . idxm - index of the matrix within the nest matrix 624d8588912SDave May . jdxm - index of the matrix within the nest matrix 625d8588912SDave May 626d8588912SDave May Output Parameter: 627d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 628d8588912SDave May 629d8588912SDave May Notes: 630d8588912SDave May 631d8588912SDave May Level: developer 632d8588912SDave May 633d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 634d8588912SDave May @*/ 635d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 636d8588912SDave May { 637699a902aSJed Brown PetscErrorCode ierr; 638d8588912SDave May 639d8588912SDave May PetscFunctionBegin; 640699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 641d8588912SDave May PetscFunctionReturn(0); 642d8588912SDave May } 643d8588912SDave May 644d8588912SDave May EXTERN_C_BEGIN 645d8588912SDave May #undef __FUNCT__ 646d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 647d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 648d8588912SDave May { 649d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 650d8588912SDave May PetscFunctionBegin; 651d8588912SDave May if (M) { *M = bA->nr; } 652d8588912SDave May if (N) { *N = bA->nc; } 653d8588912SDave May if (mat) { *mat = bA->m; } 654d8588912SDave May PetscFunctionReturn(0); 655d8588912SDave May } 656d8588912SDave May EXTERN_C_END 657d8588912SDave May 658d8588912SDave May #undef __FUNCT__ 659d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 660d8588912SDave May /*@C 661d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 662d8588912SDave May 663d8588912SDave May Not collective 664d8588912SDave May 665d8588912SDave May Input Parameters: 666d8588912SDave May . A - nest matri 667d8588912SDave May 668d8588912SDave May Output Parameter: 669d8588912SDave May . M - number of rows in the nest matrix 670d8588912SDave May . N - number of cols in the nest matrix 671d8588912SDave May . mat - 2d array of matrices 672d8588912SDave May 673d8588912SDave May Notes: 674d8588912SDave May 675d8588912SDave May The user should not free the array mat. 676d8588912SDave May 677d8588912SDave May Level: developer 678d8588912SDave May 679d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 680d8588912SDave May @*/ 681d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 682d8588912SDave May { 683699a902aSJed Brown PetscErrorCode ierr; 684d8588912SDave May 685d8588912SDave May PetscFunctionBegin; 686699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 687d8588912SDave May PetscFunctionReturn(0); 688d8588912SDave May } 689d8588912SDave May 690d8588912SDave May EXTERN_C_BEGIN 691d8588912SDave May #undef __FUNCT__ 692d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 693d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 694d8588912SDave May { 695d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 696d8588912SDave May 697d8588912SDave May PetscFunctionBegin; 698d8588912SDave May if (M) { *M = bA->nr; } 699d8588912SDave May if (N) { *N = bA->nc; } 700d8588912SDave May PetscFunctionReturn(0); 701d8588912SDave May } 702d8588912SDave May EXTERN_C_END 703d8588912SDave May 704d8588912SDave May #undef __FUNCT__ 705d8588912SDave May #define __FUNCT__ "MatNestGetSize" 706d8588912SDave May /*@C 707d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 708d8588912SDave May 709d8588912SDave May Not collective 710d8588912SDave May 711d8588912SDave May Input Parameters: 712d8588912SDave May . A - nest matrix 713d8588912SDave May 714d8588912SDave May Output Parameter: 715d8588912SDave May . M - number of rows in the nested mat 716d8588912SDave May . N - number of cols in the nested mat 717d8588912SDave May 718d8588912SDave May Notes: 719d8588912SDave May 720d8588912SDave May Level: developer 721d8588912SDave May 722d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 723d8588912SDave May @*/ 724d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 725d8588912SDave May { 726699a902aSJed Brown PetscErrorCode ierr; 727d8588912SDave May 728d8588912SDave May PetscFunctionBegin; 729699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 730d8588912SDave May PetscFunctionReturn(0); 731d8588912SDave May } 732d8588912SDave May 733207556f9SJed Brown EXTERN_C_BEGIN 734207556f9SJed Brown #undef __FUNCT__ 735207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 736207556f9SJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatNestSetVecType_Nest(Mat A,const VecType vtype) 737207556f9SJed Brown { 738207556f9SJed Brown PetscErrorCode ierr; 739207556f9SJed Brown PetscBool flg; 740207556f9SJed Brown 741207556f9SJed Brown PetscFunctionBegin; 742207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 743207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 744207556f9SJed Brown A->ops->getvecs = flg ? MatGetVecs_Nest : 0; 745207556f9SJed Brown PetscFunctionReturn(0); 746207556f9SJed Brown } 747207556f9SJed Brown EXTERN_C_END 748207556f9SJed Brown 749207556f9SJed Brown #undef __FUNCT__ 750207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 751207556f9SJed Brown /*@C 752207556f9SJed Brown MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 753207556f9SJed Brown 754207556f9SJed Brown Not collective 755207556f9SJed Brown 756207556f9SJed Brown Input Parameters: 757207556f9SJed Brown + A - nest matrix 758207556f9SJed Brown - vtype - type to use for creating vectors 759207556f9SJed Brown 760207556f9SJed Brown Notes: 761207556f9SJed Brown 762207556f9SJed Brown Level: developer 763207556f9SJed Brown 764207556f9SJed Brown .seealso: MatGetVecs() 765207556f9SJed Brown @*/ 766207556f9SJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatNestSetVecType(Mat A,const VecType vtype) 767207556f9SJed Brown { 768207556f9SJed Brown PetscErrorCode ierr; 769207556f9SJed Brown 770207556f9SJed Brown PetscFunctionBegin; 771207556f9SJed Brown ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr); 772207556f9SJed Brown PetscFunctionReturn(0); 773207556f9SJed Brown } 774207556f9SJed Brown 775d8588912SDave May /* constructor */ 776d8588912SDave May #undef __FUNCT__ 777d8588912SDave May #define __FUNCT__ "MatNestSetOps_Private" 778d8588912SDave May static PetscErrorCode MatNestSetOps_Private(struct _MatOps* ops) 779d8588912SDave May { 780d8588912SDave May PetscFunctionBegin; 781d8588912SDave May /* 0*/ 782d8588912SDave May ops->setvalues = 0; 783d8588912SDave May ops->getrow = 0; 784d8588912SDave May ops->restorerow = 0; 785d8588912SDave May ops->mult = MatMult_Nest; 786d8588912SDave May ops->multadd = 0; 787d8588912SDave May /* 5*/ 788d8588912SDave May ops->multtranspose = MatMultTranspose_Nest; 789d8588912SDave May ops->multtransposeadd = 0; 790d8588912SDave May ops->solve = 0; 791d8588912SDave May ops->solveadd = 0; 792d8588912SDave May ops->solvetranspose = 0; 793d8588912SDave May /*10*/ 794d8588912SDave May ops->solvetransposeadd = 0; 795d8588912SDave May ops->lufactor = 0; 796d8588912SDave May ops->choleskyfactor = 0; 797d8588912SDave May ops->sor = 0; 798d8588912SDave May ops->transpose = 0; 799d8588912SDave May /*15*/ 800d8588912SDave May ops->getinfo = 0; 801d8588912SDave May ops->equal = 0; 802d8588912SDave May ops->getdiagonal = 0; 803d8588912SDave May ops->diagonalscale = 0; 804d8588912SDave May ops->norm = 0; 805d8588912SDave May /*20*/ 806d8588912SDave May ops->assemblybegin = MatAssemblyBegin_Nest; 807d8588912SDave May ops->assemblyend = MatAssemblyEnd_Nest; 808d8588912SDave May ops->setoption = 0; 809d8588912SDave May ops->zeroentries = MatZeroEntries_Nest; 810d8588912SDave May /*24*/ 811d8588912SDave May ops->zerorows = 0; 812d8588912SDave May ops->lufactorsymbolic = 0; 813d8588912SDave May ops->lufactornumeric = 0; 814d8588912SDave May ops->choleskyfactorsymbolic = 0; 815d8588912SDave May ops->choleskyfactornumeric = 0; 816d8588912SDave May /*29*/ 817d8588912SDave May ops->setuppreallocation = 0; 818d8588912SDave May ops->ilufactorsymbolic = 0; 819d8588912SDave May ops->iccfactorsymbolic = 0; 820d8588912SDave May ops->getarray = 0; 821d8588912SDave May ops->restorearray = 0; 822d8588912SDave May /*34*/ 823d8588912SDave May ops->duplicate = MatDuplicate_Nest; 824d8588912SDave May ops->forwardsolve = 0; 825d8588912SDave May ops->backwardsolve = 0; 826d8588912SDave May ops->ilufactor = 0; 827d8588912SDave May ops->iccfactor = 0; 828d8588912SDave May /*39*/ 829d8588912SDave May ops->axpy = 0; 830d8588912SDave May ops->getsubmatrices = 0; 831d8588912SDave May ops->increaseoverlap = 0; 832d8588912SDave May ops->getvalues = 0; 833d8588912SDave May ops->copy = 0; 834d8588912SDave May /*44*/ 835d8588912SDave May ops->getrowmax = 0; 836d8588912SDave May ops->scale = 0; 837d8588912SDave May ops->shift = 0; 838d8588912SDave May ops->diagonalset = 0; 839d8588912SDave May ops->zerorowscolumns = 0; 840d8588912SDave May /*49*/ 841d8588912SDave May ops->setblocksize = 0; 842d8588912SDave May ops->getrowij = 0; 843d8588912SDave May ops->restorerowij = 0; 844d8588912SDave May ops->getcolumnij = 0; 845d8588912SDave May ops->restorecolumnij = 0; 846d8588912SDave May /*54*/ 847d8588912SDave May ops->fdcoloringcreate = 0; 848d8588912SDave May ops->coloringpatch = 0; 849d8588912SDave May ops->setunfactored = 0; 850d8588912SDave May ops->permute = 0; 851d8588912SDave May ops->setvaluesblocked = 0; 852d8588912SDave May /*59*/ 853f349c1fdSJed Brown ops->getsubmatrix = MatGetSubMatrix_Nest; 854d8588912SDave May ops->destroy = MatDestroy_Nest; 855d8588912SDave May ops->view = MatView_Nest; 856d8588912SDave May ops->convertfrom = 0; 857d8588912SDave May ops->usescaledform = 0; 858d8588912SDave May /*64*/ 859d8588912SDave May ops->scalesystem = 0; 860d8588912SDave May ops->unscalesystem = 0; 861d8588912SDave May ops->setlocaltoglobalmapping = 0; 862d8588912SDave May ops->setvalueslocal = 0; 863d8588912SDave May ops->zerorowslocal = 0; 864d8588912SDave May /*69*/ 865d8588912SDave May ops->getrowmaxabs = 0; 866d8588912SDave May ops->getrowminabs = 0; 867d8588912SDave May ops->convert = 0;/*MatConvert_Nest;*/ 868d8588912SDave May ops->setcoloring = 0; 869d8588912SDave May ops->setvaluesadic = 0; 870d8588912SDave May /* 74 */ 871d8588912SDave May ops->setvaluesadifor = 0; 872d8588912SDave May ops->fdcoloringapply = 0; 873d8588912SDave May ops->setfromoptions = 0; 874d8588912SDave May ops->multconstrained = 0; 875d8588912SDave May ops->multtransposeconstrained = 0; 876d8588912SDave May /*79*/ 877d8588912SDave May ops->permutesparsify = 0; 878d8588912SDave May ops->mults = 0; 879d8588912SDave May ops->solves = 0; 880d8588912SDave May ops->getinertia = 0; 881d8588912SDave May ops->load = 0; 882d8588912SDave May /*84*/ 883d8588912SDave May ops->issymmetric = 0; 884d8588912SDave May ops->ishermitian = 0; 885d8588912SDave May ops->isstructurallysymmetric = 0; 886d8588912SDave May ops->dummy4 = 0; 887207556f9SJed Brown ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 888d8588912SDave May /*89*/ 889d8588912SDave May ops->matmult = 0;/*MatMatMult_Nest;*/ 890d8588912SDave May ops->matmultsymbolic = 0; 891d8588912SDave May ops->matmultnumeric = 0; 892d8588912SDave May ops->ptap = 0; 893d8588912SDave May ops->ptapsymbolic = 0; 894d8588912SDave May /*94*/ 895d8588912SDave May ops->ptapnumeric = 0; 896d8588912SDave May ops->matmulttranspose = 0; 897d8588912SDave May ops->matmulttransposesymbolic = 0; 898d8588912SDave May ops->matmulttransposenumeric = 0; 899d8588912SDave May ops->ptapsymbolic_seqaij = 0; 900d8588912SDave May /*99*/ 901d8588912SDave May ops->ptapnumeric_seqaij = 0; 902d8588912SDave May ops->ptapsymbolic_mpiaij = 0; 903d8588912SDave May ops->ptapnumeric_mpiaij = 0; 904d8588912SDave May ops->conjugate = 0; 905d8588912SDave May ops->setsizes = 0; 906d8588912SDave May /*104*/ 907d8588912SDave May ops->setvaluesrow = 0; 908d8588912SDave May ops->realpart = 0; 909d8588912SDave May ops->imaginarypart = 0; 910d8588912SDave May ops->getrowuppertriangular = 0; 911d8588912SDave May ops->restorerowuppertriangular = 0; 912d8588912SDave May /*109*/ 913d8588912SDave May ops->matsolve = 0; 914d8588912SDave May ops->getredundantmatrix = 0; 915d8588912SDave May ops->getrowmin = 0; 916d8588912SDave May ops->getcolumnvector = 0; 917d8588912SDave May ops->missingdiagonal = 0; 918d8588912SDave May /* 114 */ 919d8588912SDave May ops->getseqnonzerostructure = 0; 920d8588912SDave May ops->create = 0; 921d8588912SDave May ops->getghosts = 0; 922d8588912SDave May ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 923d8588912SDave May ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 924d8588912SDave May /* 119 */ 925d8588912SDave May ops->multdiagonalblock = 0; 926d8588912SDave May ops->hermitiantranspose = 0; 927d8588912SDave May ops->multhermitiantranspose = 0; 928d8588912SDave May ops->multhermitiantransposeadd = 0; 929d8588912SDave May ops->getmultiprocblock = 0; 930d8588912SDave May /* 124 */ 931d8588912SDave May ops->dummy1 = 0; 932d8588912SDave May ops->dummy2 = 0; 933d8588912SDave May ops->dummy3 = 0; 934d8588912SDave May ops->dummy4 = 0; 935d8588912SDave May ops->getsubmatricesparallel = 0; 936d8588912SDave May 937d8588912SDave May PetscFunctionReturn(0); 938d8588912SDave May } 939d8588912SDave May 940d8588912SDave May #undef __FUNCT__ 941d8588912SDave May #define __FUNCT__ "MatSetUp_Nest_Private" 942841e96a3SJed Brown static PetscErrorCode MatSetUp_Nest_Private(Mat A,PetscInt nr,PetscInt nc,const Mat *sub) 943d8588912SDave May { 944d8588912SDave May Mat_Nest *ctx = (Mat_Nest*)A->data; 945d8588912SDave May PetscInt i,j; 946d8588912SDave May PetscErrorCode ierr; 947d8588912SDave May 948d8588912SDave May PetscFunctionBegin; 949d8588912SDave May if(ctx->setup_called) PetscFunctionReturn(0); 950d8588912SDave May 951d8588912SDave May ctx->nr = nr; 952d8588912SDave May ctx->nc = nc; 953d8588912SDave May 954d8588912SDave May /* Create space */ 955d8588912SDave May ierr = PetscMalloc(sizeof(Mat*)*ctx->nr,&ctx->m);CHKERRQ(ierr); 956d8588912SDave May for (i=0; i<ctx->nr; i++) { 957d8588912SDave May ierr = PetscMalloc(sizeof(Mat)*ctx->nc,&ctx->m[i]);CHKERRQ(ierr); 958d8588912SDave May } 959d8588912SDave May for (i=0; i<ctx->nr; i++) { 960d8588912SDave May for (j=0; j<ctx->nc; j++) { 961841e96a3SJed Brown ctx->m[i][j] = sub[i*nc+j]; 962841e96a3SJed Brown if (sub[i*nc+j]) { 963841e96a3SJed Brown ierr = PetscObjectReference((PetscObject)sub[i*nc+j]);CHKERRQ(ierr); 964d8588912SDave May } 965d8588912SDave May } 966d8588912SDave May } 967d8588912SDave May 968f349c1fdSJed Brown ierr = PetscMalloc(sizeof(IS)*ctx->nr,&ctx->isglobal.row);CHKERRQ(ierr); 969f349c1fdSJed Brown ierr = PetscMalloc(sizeof(IS)*ctx->nc,&ctx->isglobal.col);CHKERRQ(ierr); 970d8588912SDave May 971d8588912SDave May ierr = PetscMalloc(sizeof(PetscInt)*ctx->nr,&ctx->row_len);CHKERRQ(ierr); 972d8588912SDave May ierr = PetscMalloc(sizeof(PetscInt)*ctx->nc,&ctx->col_len);CHKERRQ(ierr); 973d8588912SDave May for (i=0;i<ctx->nr;i++) ctx->row_len[i]=-1; 974d8588912SDave May for (j=0;j<ctx->nc;j++) ctx->col_len[j]=-1; 975d8588912SDave May 976d8588912SDave May ctx->setup_called = PETSC_TRUE; 977d8588912SDave May 978d8588912SDave May PetscFunctionReturn(0); 979d8588912SDave May } 980d8588912SDave May 981d8588912SDave May 982d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 983d8588912SDave May /* 984d8588912SDave May nprocessors = NP 985d8588912SDave May Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) ) 986d8588912SDave May proc 0: => (g_0,h_0,) 987d8588912SDave May proc 1: => (g_1,h_1,) 988d8588912SDave May ... 989d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 990d8588912SDave May 991d8588912SDave May proc 0: proc 1: proc nprocs-1: 992d8588912SDave May is[0] = ( 0,1,2,...,nlocal(g_0)-1 ) ( 0,1,...,nlocal(g_1)-1 ) ( 0,1,...,nlocal(g_NP-1) ) 993d8588912SDave May 994d8588912SDave May proc 0: 995d8588912SDave May is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 ) 996d8588912SDave May proc 1: 997d8588912SDave May is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 ) 998d8588912SDave May 999d8588912SDave May proc NP-1: 1000d8588912SDave May is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 ) 1001d8588912SDave May */ 1002d8588912SDave May #undef __FUNCT__ 1003d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1004841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1005d8588912SDave May { 1006e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 10072ae74bdbSJed Brown PetscInt i,j,offset,n,bs; 1008d8588912SDave May PetscErrorCode ierr; 10092ae74bdbSJed Brown Mat sub; 1010d8588912SDave May 1011d8588912SDave May PetscFunctionBegin; 1012d8588912SDave May if (is_row) { /* valid IS is passed in */ 1013d8588912SDave May /* refs on is[] are incremeneted */ 1014e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1015d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1016e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1017d8588912SDave May } 10182ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 10192ae74bdbSJed Brown offset = A->rmap->rstart; 1020e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1021f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 10222ae74bdbSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 1023207556f9SJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 10242ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1025e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1026e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 10272ae74bdbSJed Brown offset += n; 1028d8588912SDave May } 1029d8588912SDave May } 1030d8588912SDave May 1031d8588912SDave May if (is_col) { /* valid IS is passed in */ 1032d8588912SDave May /* refs on is[] are incremeneted */ 1033e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1034d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1035e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1036d8588912SDave May } 10372ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 10382ae74bdbSJed Brown offset = A->cmap->rstart; 1039e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1040f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 10412ae74bdbSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 1042207556f9SJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 10432ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1044e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1045e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 10462ae74bdbSJed Brown offset += n; 1047d8588912SDave May } 1048d8588912SDave May } 1049e2d7f03fSJed Brown 1050e2d7f03fSJed Brown /* Set up the local ISs */ 1051e2d7f03fSJed Brown ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1052e2d7f03fSJed Brown ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1053e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1054e2d7f03fSJed Brown IS isloc; 1055e2d7f03fSJed Brown ISLocalToGlobalMapping rmap; 1056e2d7f03fSJed Brown PetscInt nlocal,bs; 1057e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1058e2d7f03fSJed Brown ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr); 1059207556f9SJed Brown if (rmap) { 1060e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1061e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1062e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1063e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1064207556f9SJed Brown } else { 1065207556f9SJed Brown nlocal = 0; 1066207556f9SJed Brown isloc = PETSC_NULL; 1067207556f9SJed Brown } 1068e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1069e2d7f03fSJed Brown offset += nlocal; 1070e2d7f03fSJed Brown } 1071e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1072e2d7f03fSJed Brown IS isloc; 1073e2d7f03fSJed Brown ISLocalToGlobalMapping cmap; 1074e2d7f03fSJed Brown PetscInt nlocal,bs; 1075e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 1076e2d7f03fSJed Brown ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr); 1077207556f9SJed Brown if (cmap) { 1078e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1079e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1080e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1081e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1082207556f9SJed Brown } else { 1083207556f9SJed Brown nlocal = 0; 1084207556f9SJed Brown isloc = PETSC_NULL; 1085207556f9SJed Brown } 1086e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1087e2d7f03fSJed Brown offset += nlocal; 1088e2d7f03fSJed Brown } 1089d8588912SDave May PetscFunctionReturn(0); 1090d8588912SDave May } 1091d8588912SDave May 1092d8588912SDave May #undef __FUNCT__ 1093d8588912SDave May #define __FUNCT__ "MatCreateNest" 1094841e96a3SJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1095d8588912SDave May { 1096d8588912SDave May Mat A; 1097d8588912SDave May Mat_Nest *s; 10982ae74bdbSJed Brown PetscInt i,m,n,M,N; 1099d8588912SDave May PetscErrorCode ierr; 1100d8588912SDave May 1101d8588912SDave May PetscFunctionBegin; 1102841e96a3SJed Brown if (nr < 0) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 11032ae74bdbSJed Brown if (nr && is_row) { 11042ae74bdbSJed Brown PetscValidPointer(is_row,3); 11052ae74bdbSJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 11062ae74bdbSJed Brown } 1107841e96a3SJed Brown if (nc < 0) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 11082ae74bdbSJed Brown if (nc && is_row) { 11092ae74bdbSJed Brown PetscValidPointer(is_col,5); 11102ae74bdbSJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 11112ae74bdbSJed Brown } 1112841e96a3SJed Brown if (nr*nc) PetscValidPointer(a,6); 1113841e96a3SJed Brown PetscValidPointer(B,7); 1114841e96a3SJed Brown 1115d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1116d8588912SDave May 1117d8588912SDave May ierr = PetscMalloc( sizeof(Mat_Nest), &s );CHKERRQ(ierr); 1118d8588912SDave May ierr = PetscMemzero( s, sizeof(Mat_Nest) );CHKERRQ(ierr); 1119d8588912SDave May A->data = (void*)s; 1120d8588912SDave May s->setup_called = PETSC_FALSE; 1121d8588912SDave May s->nr = s->nc = -1; 1122d8588912SDave May s->m = PETSC_NULL; 1123d8588912SDave May 1124d8588912SDave May /* define the operations */ 1125d8588912SDave May ierr = MatNestSetOps_Private(A->ops);CHKERRQ(ierr); 1126d8588912SDave May A->spptr = 0; 1127d8588912SDave May A->same_nonzero = PETSC_FALSE; 1128d8588912SDave May A->assembled = PETSC_FALSE; 1129d8588912SDave May 1130d8588912SDave May ierr = PetscObjectChangeTypeName((PetscObject) A, MATNEST );CHKERRQ(ierr); 1131d8588912SDave May ierr = MatSetUp_Nest_Private(A,nr,nc,a);CHKERRQ(ierr); 1132d8588912SDave May 1133d8588912SDave May m = n = 0; 1134d8588912SDave May M = N = 0; 1135d8588912SDave May ierr = MatNestGetSize_Recursive(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr); 1136d8588912SDave May ierr = MatNestGetSize_Recursive(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr); 1137d8588912SDave May 1138d8588912SDave May ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1139d8588912SDave May ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1140d8588912SDave May ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 1141d8588912SDave May ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1142d8588912SDave May ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1143d8588912SDave May ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 1144d8588912SDave May 1145d8588912SDave May ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1146d8588912SDave May ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1147d8588912SDave May 1148841e96a3SJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1149207556f9SJed Brown ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 1150d8588912SDave May 1151d8588912SDave May /* expose Nest api's */ 1152d8588912SDave May ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C","MatNestGetSubMat_Nest",MatNestGetSubMat_Nest);CHKERRQ(ierr); 1153d8588912SDave May ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C","MatNestGetSubMats_Nest",MatNestGetSubMats_Nest);CHKERRQ(ierr); 1154d8588912SDave May ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C","MatNestGetSize_Nest",MatNestGetSize_Nest);CHKERRQ(ierr); 1155207556f9SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C","MatNestSetVecType_Nest",MatNestSetVecType_Nest);CHKERRQ(ierr); 1156d8588912SDave May 1157d8588912SDave May *B = A; 1158d8588912SDave May PetscFunctionReturn(0); 1159d8588912SDave May } 1160