1d8588912SDave May 2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I "petscmat.h" I*/ 30c312b8eSJed Brown #include <petscsf.h> 4d8588912SDave May 5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]); 62a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left); 7c8883902SJed Brown 8d8588912SDave May /* private functions */ 9d8588912SDave May #undef __FUNCT__ 108188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private" 118188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N) 12d8588912SDave May { 13d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 148188e55aSJed Brown PetscInt i,j; 15d8588912SDave May PetscErrorCode ierr; 16d8588912SDave May 17d8588912SDave May PetscFunctionBegin; 188188e55aSJed Brown *m = *n = *M = *N = 0; 198188e55aSJed Brown for (i=0; i<bA->nr; i++) { /* rows */ 208188e55aSJed Brown PetscInt sm,sM; 218188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr); 228188e55aSJed Brown ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr); 238188e55aSJed Brown *m += sm; 248188e55aSJed Brown *M += sM; 25d8588912SDave May } 268188e55aSJed Brown for (j=0; j<bA->nc; j++) { /* cols */ 278188e55aSJed Brown PetscInt sn,sN; 288188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr); 298188e55aSJed Brown ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr); 308188e55aSJed Brown *n += sn; 318188e55aSJed Brown *N += sN; 32d8588912SDave May } 33d8588912SDave May PetscFunctionReturn(0); 34d8588912SDave May } 35d8588912SDave May 36d8588912SDave May /* operations */ 37d8588912SDave May #undef __FUNCT__ 38d8588912SDave May #define __FUNCT__ "MatMult_Nest" 39207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 40d8588912SDave May { 41d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 42207556f9SJed Brown Vec *bx = bA->right,*by = bA->left; 43207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 44d8588912SDave May PetscErrorCode ierr; 45d8588912SDave May 46d8588912SDave May PetscFunctionBegin; 47207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 48207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 49207556f9SJed Brown for (i=0; i<nr; i++) { 50d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 51207556f9SJed Brown for (j=0; j<nc; j++) { 52207556f9SJed Brown if (!bA->m[i][j]) continue; 53d8588912SDave May /* y[i] <- y[i] + A[i][j] * x[j] */ 54d8588912SDave May ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 55d8588912SDave May } 56d8588912SDave May } 57207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 58207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 59d8588912SDave May PetscFunctionReturn(0); 60d8588912SDave May } 61d8588912SDave May 62d8588912SDave May #undef __FUNCT__ 639194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest" 649194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z) 659194d70fSJed Brown { 669194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 679194d70fSJed Brown Vec *bx = bA->right,*bz = bA->left; 689194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 699194d70fSJed Brown PetscErrorCode ierr; 709194d70fSJed Brown 719194d70fSJed Brown PetscFunctionBegin; 729194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 739194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 749194d70fSJed Brown for (i=0; i<nr; i++) { 759194d70fSJed Brown if (y != z) { 769194d70fSJed Brown Vec by; 779194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 789194d70fSJed Brown ierr = VecCopy(by,bz[i]);CHKERRQ(ierr); 79336d21e7SJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 809194d70fSJed Brown } 819194d70fSJed Brown for (j=0; j<nc; j++) { 829194d70fSJed Brown if (!bA->m[i][j]) continue; 839194d70fSJed Brown /* y[i] <- y[i] + A[i][j] * x[j] */ 849194d70fSJed Brown ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr); 859194d70fSJed Brown } 869194d70fSJed Brown } 879194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 889194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 899194d70fSJed Brown PetscFunctionReturn(0); 909194d70fSJed Brown } 919194d70fSJed Brown 929194d70fSJed Brown #undef __FUNCT__ 93d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest" 94207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 95d8588912SDave May { 96d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 97207556f9SJed Brown Vec *bx = bA->left,*by = bA->right; 98207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 99d8588912SDave May PetscErrorCode ierr; 100d8588912SDave May 101d8588912SDave May PetscFunctionBegin; 102609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 103609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 104207556f9SJed Brown for (j=0; j<nc; j++) { 105609e31cbSJed Brown ierr = VecZeroEntries(by[j]);CHKERRQ(ierr); 106609e31cbSJed Brown for (i=0; i<nr; i++) { 1076c75ac25SJed Brown if (!bA->m[i][j]) continue; 108609e31cbSJed Brown /* y[j] <- y[j] + (A[i][j])^T * x[i] */ 109609e31cbSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr); 110d8588912SDave May } 111d8588912SDave May } 112609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 113609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 114d8588912SDave May PetscFunctionReturn(0); 115d8588912SDave May } 116d8588912SDave May 117d8588912SDave May #undef __FUNCT__ 1189194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest" 1199194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z) 1209194d70fSJed Brown { 1219194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 1229194d70fSJed Brown Vec *bx = bA->left,*bz = bA->right; 1239194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 1249194d70fSJed Brown PetscErrorCode ierr; 1259194d70fSJed Brown 1269194d70fSJed Brown PetscFunctionBegin; 1279194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 1289194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 1299194d70fSJed Brown for (j=0; j<nc; j++) { 1309194d70fSJed Brown if (y != z) { 1319194d70fSJed Brown Vec by; 1329194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 1339194d70fSJed Brown ierr = VecCopy(by,bz[j]);CHKERRQ(ierr); 1349194d70fSJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 1359194d70fSJed Brown } 1369194d70fSJed Brown for (i=0; i<nr; i++) { 1376c75ac25SJed Brown if (!bA->m[i][j]) continue; 1389194d70fSJed Brown /* z[j] <- y[j] + (A[i][j])^T * x[i] */ 1399194d70fSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr); 1409194d70fSJed Brown } 1419194d70fSJed Brown } 1429194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 1439194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 1449194d70fSJed Brown PetscFunctionReturn(0); 1459194d70fSJed Brown } 1469194d70fSJed Brown 1479194d70fSJed Brown #undef __FUNCT__ 148e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList" 149e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 150e2d7f03fSJed Brown { 151e2d7f03fSJed Brown PetscErrorCode ierr; 152e2d7f03fSJed Brown IS *lst = *list; 153e2d7f03fSJed Brown PetscInt i; 154e2d7f03fSJed Brown 155e2d7f03fSJed Brown PetscFunctionBegin; 156e2d7f03fSJed Brown if (!lst) PetscFunctionReturn(0); 1576bf464f9SBarry Smith for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);} 158e2d7f03fSJed Brown ierr = PetscFree(lst);CHKERRQ(ierr); 1590298fd71SBarry Smith *list = NULL; 160e2d7f03fSJed Brown PetscFunctionReturn(0); 161e2d7f03fSJed Brown } 162e2d7f03fSJed Brown 163e2d7f03fSJed Brown #undef __FUNCT__ 164d8588912SDave May #define __FUNCT__ "MatDestroy_Nest" 165207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A) 166d8588912SDave May { 167d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 168d8588912SDave May PetscInt i,j; 169d8588912SDave May PetscErrorCode ierr; 170d8588912SDave May 171d8588912SDave May PetscFunctionBegin; 172d8588912SDave May /* release the matrices and the place holders */ 173e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 174e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 175e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 176e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 177d8588912SDave May 178d8588912SDave May ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 179d8588912SDave May ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 180d8588912SDave May 181207556f9SJed Brown ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 182207556f9SJed Brown 183d8588912SDave May /* release the matrices and the place holders */ 184d8588912SDave May if (vs->m) { 185d8588912SDave May for (i=0; i<vs->nr; i++) { 186d8588912SDave May for (j=0; j<vs->nc; j++) { 1876bf464f9SBarry Smith ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr); 188d8588912SDave May } 189d8588912SDave May ierr = PetscFree(vs->m[i]);CHKERRQ(ierr); 190d8588912SDave May } 191d8588912SDave May ierr = PetscFree(vs->m);CHKERRQ(ierr); 192d8588912SDave May } 193bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 194d8588912SDave May 195bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr); 196bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr); 197bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr); 198bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr); 199bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr); 200bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr); 201bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr); 202bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr); 203d8588912SDave May PetscFunctionReturn(0); 204d8588912SDave May } 205d8588912SDave May 206d8588912SDave May #undef __FUNCT__ 207d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest" 208207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 209d8588912SDave May { 210d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 211d8588912SDave May PetscInt i,j; 212d8588912SDave May PetscErrorCode ierr; 213d8588912SDave May 214d8588912SDave May PetscFunctionBegin; 215d8588912SDave May for (i=0; i<vs->nr; i++) { 216d8588912SDave May for (j=0; j<vs->nc; j++) { 217e7c19651SJed Brown if (vs->m[i][j]) { 218e7c19651SJed Brown ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); 219e7c19651SJed Brown if (!vs->splitassembly) { 220e7c19651SJed Brown /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested 221e7c19651SJed Brown * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was 222e7c19651SJed Brown * already performing an assembly, but the result would by more complicated and appears to offer less 223e7c19651SJed Brown * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an 224e7c19651SJed Brown * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives". 225e7c19651SJed Brown */ 226e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 227e7c19651SJed Brown } 228e7c19651SJed Brown } 229d8588912SDave May } 230d8588912SDave May } 231d8588912SDave May PetscFunctionReturn(0); 232d8588912SDave May } 233d8588912SDave May 234d8588912SDave May #undef __FUNCT__ 235d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 236207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 237d8588912SDave May { 238d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 239d8588912SDave May PetscInt i,j; 240d8588912SDave May PetscErrorCode ierr; 241d8588912SDave May 242d8588912SDave May PetscFunctionBegin; 243d8588912SDave May for (i=0; i<vs->nr; i++) { 244d8588912SDave May for (j=0; j<vs->nc; j++) { 245e7c19651SJed Brown if (vs->m[i][j]) { 246e7c19651SJed Brown if (vs->splitassembly) { 247e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 248e7c19651SJed Brown } 249e7c19651SJed Brown } 250d8588912SDave May } 251d8588912SDave May } 252d8588912SDave May PetscFunctionReturn(0); 253d8588912SDave May } 254d8588912SDave May 255d8588912SDave May #undef __FUNCT__ 256f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 257f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 258d8588912SDave May { 259207556f9SJed Brown PetscErrorCode ierr; 260f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 261f349c1fdSJed Brown PetscInt j; 262f349c1fdSJed Brown Mat sub; 263d8588912SDave May 264d8588912SDave May PetscFunctionBegin; 2650298fd71SBarry Smith sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */ 266f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 2674994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 268f349c1fdSJed Brown *B = sub; 269f349c1fdSJed Brown PetscFunctionReturn(0); 270d8588912SDave May } 271d8588912SDave May 272f349c1fdSJed Brown #undef __FUNCT__ 273f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 274f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 275f349c1fdSJed Brown { 276207556f9SJed Brown PetscErrorCode ierr; 277f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 278f349c1fdSJed Brown PetscInt i; 279f349c1fdSJed Brown Mat sub; 280f349c1fdSJed Brown 281f349c1fdSJed Brown PetscFunctionBegin; 2820298fd71SBarry Smith sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */ 283f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 2844994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 285f349c1fdSJed Brown *B = sub; 286f349c1fdSJed Brown PetscFunctionReturn(0); 287d8588912SDave May } 288d8588912SDave May 289f349c1fdSJed Brown #undef __FUNCT__ 290f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 291f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 292f349c1fdSJed Brown { 293f349c1fdSJed Brown PetscErrorCode ierr; 294f349c1fdSJed Brown PetscInt i; 295f349c1fdSJed Brown PetscBool flg; 296f349c1fdSJed Brown 297f349c1fdSJed Brown PetscFunctionBegin; 298f349c1fdSJed Brown PetscValidPointer(list,3); 299f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 300f349c1fdSJed Brown PetscValidIntPointer(found,5); 301f349c1fdSJed Brown *found = -1; 302f349c1fdSJed Brown for (i=0; i<n; i++) { 303207556f9SJed Brown if (!list[i]) continue; 304f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 305f349c1fdSJed Brown if (flg) { 306f349c1fdSJed Brown *found = i; 307f349c1fdSJed Brown PetscFunctionReturn(0); 308f349c1fdSJed Brown } 309f349c1fdSJed Brown } 310ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set"); 311f349c1fdSJed Brown PetscFunctionReturn(0); 312f349c1fdSJed Brown } 313f349c1fdSJed Brown 314f349c1fdSJed Brown #undef __FUNCT__ 3158188e55aSJed Brown #define __FUNCT__ "MatNestGetRow" 3168188e55aSJed Brown /* Get a block row as a new MatNest */ 3178188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 3188188e55aSJed Brown { 3198188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3208188e55aSJed Brown char keyname[256]; 3218188e55aSJed Brown PetscErrorCode ierr; 3228188e55aSJed Brown 3238188e55aSJed Brown PetscFunctionBegin; 3240298fd71SBarry Smith *B = NULL; 3258caf3d72SBarry Smith ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr); 3268188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 3278188e55aSJed Brown if (*B) PetscFunctionReturn(0); 3288188e55aSJed Brown 329ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 33026fbe8dcSKarl Rupp 3318188e55aSJed Brown (*B)->assembled = A->assembled; 33226fbe8dcSKarl Rupp 3338188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 3348188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 3358188e55aSJed Brown PetscFunctionReturn(0); 3368188e55aSJed Brown } 3378188e55aSJed Brown 3388188e55aSJed Brown #undef __FUNCT__ 339f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 340f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 341f349c1fdSJed Brown { 342f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3438188e55aSJed Brown PetscErrorCode ierr; 3446b3a5b13SJed Brown PetscInt row,col; 345e072481dSJed Brown PetscBool same,isFullCol,isFullColGlobal; 346f349c1fdSJed Brown 347f349c1fdSJed Brown PetscFunctionBegin; 3488188e55aSJed Brown /* Check if full column space. This is a hack */ 3498188e55aSJed Brown isFullCol = PETSC_FALSE; 350251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 3518188e55aSJed Brown if (same) { 35277019fcaSJed Brown PetscInt n,first,step,i,an,am,afirst,astep; 3538188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 3548188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 35577019fcaSJed Brown isFullCol = PETSC_TRUE; 35605ce4453SJed Brown for (i=0,an=A->cmap->rstart; i<vs->nc; i++) { 35777019fcaSJed Brown ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr); 35877019fcaSJed Brown ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr); 35977019fcaSJed Brown if (afirst != an || astep != step) isFullCol = PETSC_FALSE; 36077019fcaSJed Brown an += am; 36177019fcaSJed Brown } 36205ce4453SJed Brown if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE; 3638188e55aSJed Brown } 364b2566f29SBarry Smith ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr); 3658188e55aSJed Brown 366*427230ceSLisandro Dalcin if (isFullColGlobal && vs->nc > 1) { 3678188e55aSJed Brown PetscInt row; 3688188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 3698188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 3708188e55aSJed Brown } else { 371f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 372f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 373b6480e04SStefano Zampini if (!vs->m[row][col]) { 374b6480e04SStefano Zampini PetscInt lr,lc; 375b6480e04SStefano Zampini 376b6480e04SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr); 377b6480e04SStefano Zampini ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr); 378b6480e04SStefano Zampini ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr); 379b6480e04SStefano Zampini ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 380b6480e04SStefano Zampini ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr); 381b6480e04SStefano Zampini ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 382b6480e04SStefano Zampini ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 383b6480e04SStefano Zampini } 384f349c1fdSJed Brown *B = vs->m[row][col]; 3858188e55aSJed Brown } 386f349c1fdSJed Brown PetscFunctionReturn(0); 387f349c1fdSJed Brown } 388f349c1fdSJed Brown 389f349c1fdSJed Brown #undef __FUNCT__ 390f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 391f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 392f349c1fdSJed Brown { 393f349c1fdSJed Brown PetscErrorCode ierr; 394f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 395f349c1fdSJed Brown Mat sub; 396f349c1fdSJed Brown 397f349c1fdSJed Brown PetscFunctionBegin; 398f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 399f349c1fdSJed Brown switch (reuse) { 400f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 4017874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 402f349c1fdSJed Brown *B = sub; 403f349c1fdSJed Brown break; 404f349c1fdSJed Brown case MAT_REUSE_MATRIX: 405ce94432eSBarry Smith if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 406f349c1fdSJed Brown break; 407f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 408f349c1fdSJed Brown break; 409511c6705SHong Zhang case MAT_INPLACE_MATRIX: /* Nothing to do */ 410511c6705SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet"); 411511c6705SHong Zhang break; 412f349c1fdSJed Brown } 413f349c1fdSJed Brown PetscFunctionReturn(0); 414f349c1fdSJed Brown } 415f349c1fdSJed Brown 416f349c1fdSJed Brown #undef __FUNCT__ 417f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 418f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 419f349c1fdSJed Brown { 420f349c1fdSJed Brown PetscErrorCode ierr; 421f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 422f349c1fdSJed Brown Mat sub; 423f349c1fdSJed Brown 424f349c1fdSJed Brown PetscFunctionBegin; 425f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 426f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 427f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 428f349c1fdSJed Brown *B = sub; 429d8588912SDave May PetscFunctionReturn(0); 430d8588912SDave May } 431d8588912SDave May 432d8588912SDave May #undef __FUNCT__ 433d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 434207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 435d8588912SDave May { 436d8588912SDave May PetscErrorCode ierr; 437f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 438f349c1fdSJed Brown Mat sub; 439d8588912SDave May 440d8588912SDave May PetscFunctionBegin; 441f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 442ce94432eSBarry Smith if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 443f349c1fdSJed Brown if (sub) { 444ce94432eSBarry Smith if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4456bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 446d8588912SDave May } 447d8588912SDave May PetscFunctionReturn(0); 448d8588912SDave May } 449d8588912SDave May 450d8588912SDave May #undef __FUNCT__ 4517874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 4527874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 4537874fa86SDave May { 4547874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4557874fa86SDave May PetscInt i; 4567874fa86SDave May PetscErrorCode ierr; 4577874fa86SDave May 4587874fa86SDave May PetscFunctionBegin; 4597874fa86SDave May for (i=0; i<bA->nr; i++) { 460429bac76SJed Brown Vec bv; 461429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4627874fa86SDave May if (bA->m[i][i]) { 463429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 4647874fa86SDave May } else { 4655159a857SMatthew G. Knepley ierr = VecSet(bv,0.0);CHKERRQ(ierr); 4667874fa86SDave May } 467429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4687874fa86SDave May } 4697874fa86SDave May PetscFunctionReturn(0); 4707874fa86SDave May } 4717874fa86SDave May 4727874fa86SDave May #undef __FUNCT__ 4737874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 4747874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 4757874fa86SDave May { 4767874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 477429bac76SJed Brown Vec bl,*br; 4787874fa86SDave May PetscInt i,j; 4797874fa86SDave May PetscErrorCode ierr; 4807874fa86SDave May 4817874fa86SDave May PetscFunctionBegin; 4823f800ebeSJed Brown ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr); 4832e6472ebSElliott Sales de Andrade if (r) { 484429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 4852e6472ebSElliott Sales de Andrade } 4862e6472ebSElliott Sales de Andrade bl = NULL; 4877874fa86SDave May for (i=0; i<bA->nr; i++) { 4882e6472ebSElliott Sales de Andrade if (l) { 489429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4902e6472ebSElliott Sales de Andrade } 4917874fa86SDave May for (j=0; j<bA->nc; j++) { 4927874fa86SDave May if (bA->m[i][j]) { 493429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 4947874fa86SDave May } 4957874fa86SDave May } 4962e6472ebSElliott Sales de Andrade if (l) { 497a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4987874fa86SDave May } 4992e6472ebSElliott Sales de Andrade } 5002e6472ebSElliott Sales de Andrade if (r) { 501429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 5022e6472ebSElliott Sales de Andrade } 503429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 5047874fa86SDave May PetscFunctionReturn(0); 5057874fa86SDave May } 5067874fa86SDave May 5077874fa86SDave May #undef __FUNCT__ 508a061e289SJed Brown #define __FUNCT__ "MatScale_Nest" 509a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 510a061e289SJed Brown { 511a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 512a061e289SJed Brown PetscInt i,j; 513a061e289SJed Brown PetscErrorCode ierr; 514a061e289SJed Brown 515a061e289SJed Brown PetscFunctionBegin; 516a061e289SJed Brown for (i=0; i<bA->nr; i++) { 517a061e289SJed Brown for (j=0; j<bA->nc; j++) { 518a061e289SJed Brown if (bA->m[i][j]) { 519a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 520a061e289SJed Brown } 521a061e289SJed Brown } 522a061e289SJed Brown } 523a061e289SJed Brown PetscFunctionReturn(0); 524a061e289SJed Brown } 525a061e289SJed Brown 526a061e289SJed Brown #undef __FUNCT__ 527a061e289SJed Brown #define __FUNCT__ "MatShift_Nest" 528a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 529a061e289SJed Brown { 530a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 531a061e289SJed Brown PetscInt i; 532a061e289SJed Brown PetscErrorCode ierr; 533a061e289SJed Brown 534a061e289SJed Brown PetscFunctionBegin; 535a061e289SJed Brown for (i=0; i<bA->nr; i++) { 536ce94432eSBarry Smith if (!bA->m[i][i]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i); 537a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 538a061e289SJed Brown } 539a061e289SJed Brown PetscFunctionReturn(0); 540a061e289SJed Brown } 541a061e289SJed Brown 542a061e289SJed Brown #undef __FUNCT__ 5432a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest" 5442a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left) 545d8588912SDave May { 546d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 547d8588912SDave May Vec *L,*R; 548d8588912SDave May MPI_Comm comm; 549d8588912SDave May PetscInt i,j; 550d8588912SDave May PetscErrorCode ierr; 551d8588912SDave May 552d8588912SDave May PetscFunctionBegin; 553ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 554d8588912SDave May if (right) { 555d8588912SDave May /* allocate R */ 556854ce69bSBarry Smith ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr); 557d8588912SDave May /* Create the right vectors */ 558d8588912SDave May for (j=0; j<bA->nc; j++) { 559d8588912SDave May for (i=0; i<bA->nr; i++) { 560d8588912SDave May if (bA->m[i][j]) { 5612a7a6963SBarry Smith ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr); 562d8588912SDave May break; 563d8588912SDave May } 564d8588912SDave May } 5656c4ed002SBarry Smith if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 566d8588912SDave May } 567f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 568d8588912SDave May /* hand back control to the nest vector */ 569d8588912SDave May for (j=0; j<bA->nc; j++) { 5706bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 571d8588912SDave May } 572d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 573d8588912SDave May } 574d8588912SDave May 575d8588912SDave May if (left) { 576d8588912SDave May /* allocate L */ 577854ce69bSBarry Smith ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr); 578d8588912SDave May /* Create the left vectors */ 579d8588912SDave May for (i=0; i<bA->nr; i++) { 580d8588912SDave May for (j=0; j<bA->nc; j++) { 581d8588912SDave May if (bA->m[i][j]) { 5822a7a6963SBarry Smith ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr); 583d8588912SDave May break; 584d8588912SDave May } 585d8588912SDave May } 5866c4ed002SBarry Smith if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 587d8588912SDave May } 588d8588912SDave May 589f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 590d8588912SDave May for (i=0; i<bA->nr; i++) { 5916bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 592d8588912SDave May } 593d8588912SDave May 594d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 595d8588912SDave May } 596d8588912SDave May PetscFunctionReturn(0); 597d8588912SDave May } 598d8588912SDave May 599d8588912SDave May #undef __FUNCT__ 600d8588912SDave May #define __FUNCT__ "MatView_Nest" 601207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 602d8588912SDave May { 603d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 604d8588912SDave May PetscBool isascii; 605d8588912SDave May PetscInt i,j; 606d8588912SDave May PetscErrorCode ierr; 607d8588912SDave May 608d8588912SDave May PetscFunctionBegin; 609251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 610d8588912SDave May if (isascii) { 611d8588912SDave May 612d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n"); 613d8588912SDave May PetscViewerASCIIPushTab(viewer); /* push0 */ 614d8588912SDave May PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 615d8588912SDave May 616d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n"); 617d8588912SDave May for (i=0; i<bA->nr; i++) { 618d8588912SDave May for (j=0; j<bA->nc; j++) { 61919fd82e9SBarry Smith MatType type; 620270f95d7SJed Brown char name[256] = "",prefix[256] = ""; 621d8588912SDave May PetscInt NR,NC; 622d8588912SDave May PetscBool isNest = PETSC_FALSE; 623d8588912SDave May 624d8588912SDave May if (!bA->m[i][j]) { 6250298fd71SBarry Smith PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j); 626d8588912SDave May continue; 627d8588912SDave May } 628d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 629d8588912SDave May ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr); 6308caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 6318caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 632251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 633d8588912SDave May 634270f95d7SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 635d8588912SDave May 636d8588912SDave May if (isNest) { 637270f95d7SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 638d8588912SDave May ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 639270f95d7SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 640d8588912SDave May } 641d8588912SDave May } 642d8588912SDave May } 643d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 644d8588912SDave May } 645d8588912SDave May PetscFunctionReturn(0); 646d8588912SDave May } 647d8588912SDave May 648d8588912SDave May #undef __FUNCT__ 649d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 650207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 651d8588912SDave May { 652d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 653d8588912SDave May PetscInt i,j; 654d8588912SDave May PetscErrorCode ierr; 655d8588912SDave May 656d8588912SDave May PetscFunctionBegin; 657d8588912SDave May for (i=0; i<bA->nr; i++) { 658d8588912SDave May for (j=0; j<bA->nc; j++) { 659d8588912SDave May if (!bA->m[i][j]) continue; 660d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 661d8588912SDave May } 662d8588912SDave May } 663d8588912SDave May PetscFunctionReturn(0); 664d8588912SDave May } 665d8588912SDave May 666d8588912SDave May #undef __FUNCT__ 667c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest" 668c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str) 669c222c20dSDavid Ham { 670c222c20dSDavid Ham Mat_Nest *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data; 671c222c20dSDavid Ham PetscInt i,j,nr = bA->nr,nc = bA->nc; 672c222c20dSDavid Ham PetscErrorCode ierr; 673c222c20dSDavid Ham 674c222c20dSDavid Ham PetscFunctionBegin; 675c222c20dSDavid Ham if (nr != bB->nr || nc != bB->nc) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot copy a Mat_Nest of block size (%D,%D) to a Mat_Nest of block size (%D,%D)",bB->nr,bB->nc,nr,nc); 676c222c20dSDavid Ham for (i=0; i<nr; i++) { 677c222c20dSDavid Ham for (j=0; j<nc; j++) { 67846a2b97cSJed Brown if (bA->m[i][j] && bB->m[i][j]) { 679c222c20dSDavid Ham ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr); 68046a2b97cSJed Brown } else if (bA->m[i][j] || bB->m[i][j]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j); 681c222c20dSDavid Ham } 682c222c20dSDavid Ham } 683c222c20dSDavid Ham PetscFunctionReturn(0); 684c222c20dSDavid Ham } 685c222c20dSDavid Ham 686c222c20dSDavid Ham #undef __FUNCT__ 687d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 688207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 689d8588912SDave May { 690d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 691841e96a3SJed Brown Mat *b; 692841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 693d8588912SDave May PetscErrorCode ierr; 694d8588912SDave May 695d8588912SDave May PetscFunctionBegin; 696785e854fSJed Brown ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr); 697841e96a3SJed Brown for (i=0; i<nr; i++) { 698841e96a3SJed Brown for (j=0; j<nc; j++) { 699841e96a3SJed Brown if (bA->m[i][j]) { 700841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 701841e96a3SJed Brown } else { 7020298fd71SBarry Smith b[i*nc+j] = NULL; 703d8588912SDave May } 704d8588912SDave May } 705d8588912SDave May } 706ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 707841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 708841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 7096bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 710d8588912SDave May } 711d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 712d8588912SDave May 713841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 714841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 715d8588912SDave May PetscFunctionReturn(0); 716d8588912SDave May } 717d8588912SDave May 718d8588912SDave May /* nest api */ 719d8588912SDave May #undef __FUNCT__ 720d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 721d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 722d8588912SDave May { 723d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 7245fd66863SKarl Rupp 725d8588912SDave May PetscFunctionBegin; 726ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 727ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 728d8588912SDave May *mat = bA->m[idxm][jdxm]; 729d8588912SDave May PetscFunctionReturn(0); 730d8588912SDave May } 731d8588912SDave May 732d8588912SDave May #undef __FUNCT__ 733d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 7349ba0d327SJed Brown /*@ 735d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 736d8588912SDave May 737d8588912SDave May Not collective 738d8588912SDave May 739d8588912SDave May Input Parameters: 740629881c0SJed Brown + A - nest matrix 741d8588912SDave May . idxm - index of the matrix within the nest matrix 742629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 743d8588912SDave May 744d8588912SDave May Output Parameter: 745d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 746d8588912SDave May 747d8588912SDave May Level: developer 748d8588912SDave May 749d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 750d8588912SDave May @*/ 7517087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 752d8588912SDave May { 753699a902aSJed Brown PetscErrorCode ierr; 754d8588912SDave May 755d8588912SDave May PetscFunctionBegin; 756699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 757d8588912SDave May PetscFunctionReturn(0); 758d8588912SDave May } 759d8588912SDave May 760d8588912SDave May #undef __FUNCT__ 7610782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 7620782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 7630782ca92SJed Brown { 7640782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 7650782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 7660782ca92SJed Brown PetscErrorCode ierr; 7670782ca92SJed Brown 7680782ca92SJed Brown PetscFunctionBegin; 769ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 770ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 7710782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 7720782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 7730782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 7740782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 7750782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 7760782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 777ce94432eSBarry Smith if (M != Mi || N != Ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni); 778ce94432eSBarry Smith if (m != mi || n != ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni); 77926fbe8dcSKarl Rupp 7800782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 7810782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 7820782ca92SJed Brown bA->m[idxm][jdxm] = mat; 7830782ca92SJed Brown PetscFunctionReturn(0); 7840782ca92SJed Brown } 7850782ca92SJed Brown 7860782ca92SJed Brown #undef __FUNCT__ 7870782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 7889ba0d327SJed Brown /*@ 7890782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 7900782ca92SJed Brown 7910782ca92SJed Brown Logically collective on the submatrix communicator 7920782ca92SJed Brown 7930782ca92SJed Brown Input Parameters: 7940782ca92SJed Brown + A - nest matrix 7950782ca92SJed Brown . idxm - index of the matrix within the nest matrix 7960782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 7970782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 7980782ca92SJed Brown 7990782ca92SJed Brown Notes: 8000782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 8010782ca92SJed Brown 8020782ca92SJed Brown This increments the reference count of the submatrix. 8030782ca92SJed Brown 8040782ca92SJed Brown Level: developer 8050782ca92SJed Brown 8060782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 8070782ca92SJed Brown @*/ 8080782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 8090782ca92SJed Brown { 8100782ca92SJed Brown PetscErrorCode ierr; 8110782ca92SJed Brown 8120782ca92SJed Brown PetscFunctionBegin; 8130782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 8140782ca92SJed Brown PetscFunctionReturn(0); 8150782ca92SJed Brown } 8160782ca92SJed Brown 8170782ca92SJed Brown #undef __FUNCT__ 818d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 819d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 820d8588912SDave May { 821d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 8225fd66863SKarl Rupp 823d8588912SDave May PetscFunctionBegin; 82426fbe8dcSKarl Rupp if (M) *M = bA->nr; 82526fbe8dcSKarl Rupp if (N) *N = bA->nc; 82626fbe8dcSKarl Rupp if (mat) *mat = bA->m; 827d8588912SDave May PetscFunctionReturn(0); 828d8588912SDave May } 829d8588912SDave May 830d8588912SDave May #undef __FUNCT__ 831d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 832d8588912SDave May /*@C 833d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 834d8588912SDave May 835d8588912SDave May Not collective 836d8588912SDave May 837d8588912SDave May Input Parameters: 838629881c0SJed Brown . A - nest matrix 839d8588912SDave May 840d8588912SDave May Output Parameter: 841629881c0SJed Brown + M - number of rows in the nest matrix 842d8588912SDave May . N - number of cols in the nest matrix 843629881c0SJed Brown - mat - 2d array of matrices 844d8588912SDave May 845d8588912SDave May Notes: 846d8588912SDave May 847d8588912SDave May The user should not free the array mat. 848d8588912SDave May 849d8588912SDave May Level: developer 850d8588912SDave May 851d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 852d8588912SDave May @*/ 8537087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 854d8588912SDave May { 855699a902aSJed Brown PetscErrorCode ierr; 856d8588912SDave May 857d8588912SDave May PetscFunctionBegin; 858699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 859d8588912SDave May PetscFunctionReturn(0); 860d8588912SDave May } 861d8588912SDave May 862d8588912SDave May #undef __FUNCT__ 863d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 8647087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 865d8588912SDave May { 866d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 867d8588912SDave May 868d8588912SDave May PetscFunctionBegin; 86926fbe8dcSKarl Rupp if (M) *M = bA->nr; 87026fbe8dcSKarl Rupp if (N) *N = bA->nc; 871d8588912SDave May PetscFunctionReturn(0); 872d8588912SDave May } 873d8588912SDave May 874d8588912SDave May #undef __FUNCT__ 875d8588912SDave May #define __FUNCT__ "MatNestGetSize" 8769ba0d327SJed Brown /*@ 877d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 878d8588912SDave May 879d8588912SDave May Not collective 880d8588912SDave May 881d8588912SDave May Input Parameters: 882d8588912SDave May . A - nest matrix 883d8588912SDave May 884d8588912SDave May Output Parameter: 885629881c0SJed Brown + M - number of rows in the nested mat 886629881c0SJed Brown - N - number of cols in the nested mat 887d8588912SDave May 888d8588912SDave May Notes: 889d8588912SDave May 890d8588912SDave May Level: developer 891d8588912SDave May 892d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 893d8588912SDave May @*/ 8947087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 895d8588912SDave May { 896699a902aSJed Brown PetscErrorCode ierr; 897d8588912SDave May 898d8588912SDave May PetscFunctionBegin; 899699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 900d8588912SDave May PetscFunctionReturn(0); 901d8588912SDave May } 902d8588912SDave May 903900e7ff2SJed Brown #undef __FUNCT__ 904900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest" 905f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 906900e7ff2SJed Brown { 907900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 908900e7ff2SJed Brown PetscInt i; 909900e7ff2SJed Brown 910900e7ff2SJed Brown PetscFunctionBegin; 911900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 912900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 913900e7ff2SJed Brown PetscFunctionReturn(0); 914900e7ff2SJed Brown } 915900e7ff2SJed Brown 916900e7ff2SJed Brown #undef __FUNCT__ 917900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs" 9183a4d7b9aSSatish Balay /*@C 919900e7ff2SJed Brown MatNestGetISs - Returns the index sets partitioning the row and column spaces 920900e7ff2SJed Brown 921900e7ff2SJed Brown Not collective 922900e7ff2SJed Brown 923900e7ff2SJed Brown Input Parameters: 924900e7ff2SJed Brown . A - nest matrix 925900e7ff2SJed Brown 926900e7ff2SJed Brown Output Parameter: 927900e7ff2SJed Brown + rows - array of row index sets 928900e7ff2SJed Brown - cols - array of column index sets 929900e7ff2SJed Brown 930900e7ff2SJed Brown Level: advanced 931900e7ff2SJed Brown 932900e7ff2SJed Brown Notes: 933900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 934900e7ff2SJed Brown 935900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs() 936900e7ff2SJed Brown @*/ 937900e7ff2SJed Brown PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 938900e7ff2SJed Brown { 939900e7ff2SJed Brown PetscErrorCode ierr; 940900e7ff2SJed Brown 941900e7ff2SJed Brown PetscFunctionBegin; 942900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 943900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 944900e7ff2SJed Brown PetscFunctionReturn(0); 945900e7ff2SJed Brown } 946900e7ff2SJed Brown 947900e7ff2SJed Brown #undef __FUNCT__ 948900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest" 949f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 950900e7ff2SJed Brown { 951900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 952900e7ff2SJed Brown PetscInt i; 953900e7ff2SJed Brown 954900e7ff2SJed Brown PetscFunctionBegin; 955900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 956900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 957900e7ff2SJed Brown PetscFunctionReturn(0); 958900e7ff2SJed Brown } 959900e7ff2SJed Brown 960900e7ff2SJed Brown #undef __FUNCT__ 961900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs" 962900e7ff2SJed Brown /*@C 963900e7ff2SJed Brown MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 964900e7ff2SJed Brown 965900e7ff2SJed Brown Not collective 966900e7ff2SJed Brown 967900e7ff2SJed Brown Input Parameters: 968900e7ff2SJed Brown . A - nest matrix 969900e7ff2SJed Brown 970900e7ff2SJed Brown Output Parameter: 9710298fd71SBarry Smith + rows - array of row index sets (or NULL to ignore) 9720298fd71SBarry Smith - cols - array of column index sets (or NULL to ignore) 973900e7ff2SJed Brown 974900e7ff2SJed Brown Level: advanced 975900e7ff2SJed Brown 976900e7ff2SJed Brown Notes: 977900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 978900e7ff2SJed Brown 979900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs() 980900e7ff2SJed Brown @*/ 981900e7ff2SJed Brown PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 982900e7ff2SJed Brown { 983900e7ff2SJed Brown PetscErrorCode ierr; 984900e7ff2SJed Brown 985900e7ff2SJed Brown PetscFunctionBegin; 986900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 987900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 988900e7ff2SJed Brown PetscFunctionReturn(0); 989900e7ff2SJed Brown } 990900e7ff2SJed Brown 991207556f9SJed Brown #undef __FUNCT__ 992207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 99319fd82e9SBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,VecType vtype) 994207556f9SJed Brown { 995207556f9SJed Brown PetscErrorCode ierr; 996207556f9SJed Brown PetscBool flg; 997207556f9SJed Brown 998207556f9SJed Brown PetscFunctionBegin; 999207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 1000207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 10012a7a6963SBarry Smith if (flg) A->ops->getvecs = MatCreateVecs_Nest; 100212b53f24SSatish Balay else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0; 1003207556f9SJed Brown PetscFunctionReturn(0); 1004207556f9SJed Brown } 1005207556f9SJed Brown 1006207556f9SJed Brown #undef __FUNCT__ 1007207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 1008207556f9SJed Brown /*@C 10092a7a6963SBarry Smith MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs() 1010207556f9SJed Brown 1011207556f9SJed Brown Not collective 1012207556f9SJed Brown 1013207556f9SJed Brown Input Parameters: 1014207556f9SJed Brown + A - nest matrix 1015207556f9SJed Brown - vtype - type to use for creating vectors 1016207556f9SJed Brown 1017207556f9SJed Brown Notes: 1018207556f9SJed Brown 1019207556f9SJed Brown Level: developer 1020207556f9SJed Brown 10212a7a6963SBarry Smith .seealso: MatCreateVecs() 1022207556f9SJed Brown @*/ 102319fd82e9SBarry Smith PetscErrorCode MatNestSetVecType(Mat A,VecType vtype) 1024207556f9SJed Brown { 1025207556f9SJed Brown PetscErrorCode ierr; 1026207556f9SJed Brown 1027207556f9SJed Brown PetscFunctionBegin; 102819fd82e9SBarry Smith ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr); 1029207556f9SJed Brown PetscFunctionReturn(0); 1030207556f9SJed Brown } 1031207556f9SJed Brown 1032d8588912SDave May #undef __FUNCT__ 1033c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 1034c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1035d8588912SDave May { 1036c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 1037c8883902SJed Brown PetscInt i,j,m,n,M,N; 1038d8588912SDave May PetscErrorCode ierr; 1039d8588912SDave May 1040d8588912SDave May PetscFunctionBegin; 1041c8883902SJed Brown s->nr = nr; 1042c8883902SJed Brown s->nc = nc; 1043d8588912SDave May 1044c8883902SJed Brown /* Create space for submatrices */ 1045854ce69bSBarry Smith ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr); 1046c8883902SJed Brown for (i=0; i<nr; i++) { 1047854ce69bSBarry Smith ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr); 1048d8588912SDave May } 1049c8883902SJed Brown for (i=0; i<nr; i++) { 1050c8883902SJed Brown for (j=0; j<nc; j++) { 1051c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1052c8883902SJed Brown if (a[i*nc+j]) { 1053c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1054d8588912SDave May } 1055d8588912SDave May } 1056d8588912SDave May } 1057d8588912SDave May 10588188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1059d8588912SDave May 1060854ce69bSBarry Smith ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr); 1061854ce69bSBarry Smith ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr); 1062c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1063c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1064d8588912SDave May 10658188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1066d8588912SDave May 1067c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1068c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1069c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1070c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1071c8883902SJed Brown 1072c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1073c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1074c8883902SJed Brown 10751795a4d1SJed Brown ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr); 1076d8588912SDave May PetscFunctionReturn(0); 1077d8588912SDave May } 1078d8588912SDave May 1079c8883902SJed Brown #undef __FUNCT__ 1080c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 1081c8883902SJed Brown /*@ 1082c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1083c8883902SJed Brown 1084c8883902SJed Brown Collective on Mat 1085c8883902SJed Brown 1086c8883902SJed Brown Input Parameter: 1087c8883902SJed Brown + N - nested matrix 1088c8883902SJed Brown . nr - number of nested row blocks 10890298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1090c8883902SJed Brown . nc - number of nested column blocks 10910298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 10920298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1093c8883902SJed Brown 1094c8883902SJed Brown Level: advanced 1095c8883902SJed Brown 1096c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 1097c8883902SJed Brown @*/ 1098c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1099c8883902SJed Brown { 1100c8883902SJed Brown PetscErrorCode ierr; 1101c8883902SJed Brown PetscInt i; 1102c8883902SJed Brown 1103c8883902SJed Brown PetscFunctionBegin; 1104c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1105ce94432eSBarry Smith if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1106c8883902SJed Brown if (nr && is_row) { 1107c8883902SJed Brown PetscValidPointer(is_row,3); 1108c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1109c8883902SJed Brown } 1110ce94432eSBarry Smith if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 11111664e352SJed Brown if (nc && is_col) { 1112c8883902SJed Brown PetscValidPointer(is_col,5); 11139b30a8f6SBarry Smith for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1114c8883902SJed Brown } 1115c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 1116c8883902SJed 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); 1117c8883902SJed Brown PetscFunctionReturn(0); 1118c8883902SJed Brown } 1119d8588912SDave May 112077019fcaSJed Brown #undef __FUNCT__ 112177019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private" 112245b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog) 112377019fcaSJed Brown { 112477019fcaSJed Brown PetscErrorCode ierr; 112577019fcaSJed Brown PetscBool flg; 112677019fcaSJed Brown PetscInt i,j,m,mi,*ix; 112777019fcaSJed Brown 112877019fcaSJed Brown PetscFunctionBegin; 112977019fcaSJed Brown for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 113077019fcaSJed Brown if (islocal[i]) { 113177019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 113277019fcaSJed Brown flg = PETSC_TRUE; /* We found a non-trivial entry */ 113377019fcaSJed Brown } else { 113477019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 113577019fcaSJed Brown } 113677019fcaSJed Brown m += mi; 113777019fcaSJed Brown } 113877019fcaSJed Brown if (flg) { 1139785e854fSJed Brown ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr); 114077019fcaSJed Brown for (i=0,n=0; i<n; i++) { 11410298fd71SBarry Smith ISLocalToGlobalMapping smap = NULL; 114277019fcaSJed Brown VecScatter scat; 114377019fcaSJed Brown IS isreq; 114477019fcaSJed Brown Vec lvec,gvec; 11453361c9a7SJed Brown union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x; 114677019fcaSJed Brown Mat sub; 114777019fcaSJed Brown 1148ce94432eSBarry Smith if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers"); 114977019fcaSJed Brown if (colflg) { 115077019fcaSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 115177019fcaSJed Brown } else { 115277019fcaSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 115377019fcaSJed Brown } 11540298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);} 115577019fcaSJed Brown if (islocal[i]) { 115677019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 115777019fcaSJed Brown } else { 115877019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 115977019fcaSJed Brown } 116077019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = j; 116177019fcaSJed Brown if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);} 116277019fcaSJed Brown /* 116377019fcaSJed Brown Now we need to extract the monolithic global indices that correspond to the given split global indices. 116477019fcaSJed Brown In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces. 116577019fcaSJed Brown The approach here is ugly because it uses VecScatter to move indices. 116677019fcaSJed Brown */ 116777019fcaSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr); 116877019fcaSJed Brown ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr); 116977019fcaSJed Brown ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr); 11700298fd71SBarry Smith ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr); 117177019fcaSJed Brown ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 117277019fcaSJed Brown for (j=0; j<mi; j++) x[j].integer = ix[m+j]; 117377019fcaSJed Brown ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 117477019fcaSJed Brown ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 117577019fcaSJed Brown ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 117677019fcaSJed Brown ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 117777019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = x[j].integer; 117877019fcaSJed Brown ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 117977019fcaSJed Brown ierr = VecDestroy(&lvec);CHKERRQ(ierr); 118077019fcaSJed Brown ierr = VecDestroy(&gvec);CHKERRQ(ierr); 118177019fcaSJed Brown ierr = ISDestroy(&isreq);CHKERRQ(ierr); 118277019fcaSJed Brown ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 118377019fcaSJed Brown m += mi; 118477019fcaSJed Brown } 1185f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 118677019fcaSJed Brown } else { 11870298fd71SBarry Smith *ltog = NULL; 118877019fcaSJed Brown } 118977019fcaSJed Brown PetscFunctionReturn(0); 119077019fcaSJed Brown } 119177019fcaSJed Brown 119277019fcaSJed Brown 1193d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1194d8588912SDave May /* 1195d8588912SDave May nprocessors = NP 1196d8588912SDave May Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1)) 1197d8588912SDave May proc 0: => (g_0,h_0,) 1198d8588912SDave May proc 1: => (g_1,h_1,) 1199d8588912SDave May ... 1200d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1201d8588912SDave May 1202d8588912SDave May proc 0: proc 1: proc nprocs-1: 1203d8588912SDave May is[0] = (0,1,2,...,nlocal(g_0)-1) (0,1,...,nlocal(g_1)-1) (0,1,...,nlocal(g_NP-1)) 1204d8588912SDave May 1205d8588912SDave May proc 0: 1206d8588912SDave May is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1) 1207d8588912SDave May proc 1: 1208d8588912SDave May is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1) 1209d8588912SDave May 1210d8588912SDave May proc NP-1: 1211d8588912SDave May is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1) 1212d8588912SDave May */ 1213d8588912SDave May #undef __FUNCT__ 1214d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1215841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1216d8588912SDave May { 1217e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 12188188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1219d8588912SDave May PetscErrorCode ierr; 12200298fd71SBarry Smith Mat sub = NULL; 1221d8588912SDave May 1222d8588912SDave May PetscFunctionBegin; 1223854ce69bSBarry Smith ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr); 1224854ce69bSBarry Smith ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr); 1225d8588912SDave May if (is_row) { /* valid IS is passed in */ 1226d8588912SDave May /* refs on is[] are incremeneted */ 1227e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1228d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 122926fbe8dcSKarl Rupp 1230e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1231d8588912SDave May } 12322ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 12338188e55aSJed Brown nsum = 0; 12348188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 12358188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1236ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 12370298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 1238ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12398188e55aSJed Brown nsum += n; 12408188e55aSJed Brown } 1241ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 124230bc264bSJed Brown offset -= nsum; 1243e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1244f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12450298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 12462ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1247ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1248e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 12492ae74bdbSJed Brown offset += n; 1250d8588912SDave May } 1251d8588912SDave May } 1252d8588912SDave May 1253d8588912SDave May if (is_col) { /* valid IS is passed in */ 1254d8588912SDave May /* refs on is[] are incremeneted */ 1255e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1256d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 125726fbe8dcSKarl Rupp 1258e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1259d8588912SDave May } 12602ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 12612ae74bdbSJed Brown offset = A->cmap->rstart; 12628188e55aSJed Brown nsum = 0; 12638188e55aSJed Brown for (j=0; j<vs->nc; j++) { 12648188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1265ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 12660298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 1267ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12688188e55aSJed Brown nsum += n; 12698188e55aSJed Brown } 1270ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 127130bc264bSJed Brown offset -= nsum; 1272e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1273f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 12740298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 12752ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1276ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1277e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 12782ae74bdbSJed Brown offset += n; 1279d8588912SDave May } 1280d8588912SDave May } 1281e2d7f03fSJed Brown 1282e2d7f03fSJed Brown /* Set up the local ISs */ 1283785e854fSJed Brown ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 1284785e854fSJed Brown ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 1285e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1286e2d7f03fSJed Brown IS isloc; 12870298fd71SBarry Smith ISLocalToGlobalMapping rmap = NULL; 1288e2d7f03fSJed Brown PetscInt nlocal,bs; 1289e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12900298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);} 1291207556f9SJed Brown if (rmap) { 1292e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1293e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1294e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1295e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1296207556f9SJed Brown } else { 1297207556f9SJed Brown nlocal = 0; 12980298fd71SBarry Smith isloc = NULL; 1299207556f9SJed Brown } 1300e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1301e2d7f03fSJed Brown offset += nlocal; 1302e2d7f03fSJed Brown } 13038188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1304e2d7f03fSJed Brown IS isloc; 13050298fd71SBarry Smith ISLocalToGlobalMapping cmap = NULL; 1306e2d7f03fSJed Brown PetscInt nlocal,bs; 1307e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 13080298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);} 1309207556f9SJed Brown if (cmap) { 1310e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1311e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1312e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1313e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1314207556f9SJed Brown } else { 1315207556f9SJed Brown nlocal = 0; 13160298fd71SBarry Smith isloc = NULL; 1317207556f9SJed Brown } 1318e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1319e2d7f03fSJed Brown offset += nlocal; 1320e2d7f03fSJed Brown } 13210189643fSJed Brown 132277019fcaSJed Brown /* Set up the aggregate ISLocalToGlobalMapping */ 132377019fcaSJed Brown { 132445b6f7e9SBarry Smith ISLocalToGlobalMapping rmap,cmap; 132545b6f7e9SBarry Smith ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr); 132645b6f7e9SBarry Smith ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr); 132777019fcaSJed Brown if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 132877019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 132977019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 133077019fcaSJed Brown } 133177019fcaSJed Brown 13320189643fSJed Brown #if defined(PETSC_USE_DEBUG) 13330189643fSJed Brown for (i=0; i<vs->nr; i++) { 13340189643fSJed Brown for (j=0; j<vs->nc; j++) { 13350189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 13360189643fSJed Brown Mat B = vs->m[i][j]; 13370189643fSJed Brown if (!B) continue; 13380189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 13390189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 13400189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 13410189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 13420189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 13430189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 1344ce94432eSBarry Smith if (M != Mi || N != Ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni); 1345ce94432eSBarry Smith if (m != mi || n != ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni); 13460189643fSJed Brown } 13470189643fSJed Brown } 13480189643fSJed Brown #endif 1349a061e289SJed Brown 1350a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1351a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1352a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1353a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1354a061e289SJed Brown } 1355a061e289SJed Brown } 1356a061e289SJed Brown A->assembled = PETSC_TRUE; 1357d8588912SDave May PetscFunctionReturn(0); 1358d8588912SDave May } 1359d8588912SDave May 1360d8588912SDave May #undef __FUNCT__ 1361d8588912SDave May #define __FUNCT__ "MatCreateNest" 136245c38901SJed Brown /*@C 1363659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1364659c6bb0SJed Brown 1365659c6bb0SJed Brown Collective on Mat 1366659c6bb0SJed Brown 1367659c6bb0SJed Brown Input Parameter: 1368659c6bb0SJed Brown + comm - Communicator for the new Mat 1369659c6bb0SJed Brown . nr - number of nested row blocks 13700298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1371659c6bb0SJed Brown . nc - number of nested column blocks 13720298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 13730298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1374659c6bb0SJed Brown 1375659c6bb0SJed Brown Output Parameter: 1376659c6bb0SJed Brown . B - new matrix 1377659c6bb0SJed Brown 1378659c6bb0SJed Brown Level: advanced 1379659c6bb0SJed Brown 1380950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST 1381659c6bb0SJed Brown @*/ 13827087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1383d8588912SDave May { 1384d8588912SDave May Mat A; 1385d8588912SDave May PetscErrorCode ierr; 1386d8588912SDave May 1387d8588912SDave May PetscFunctionBegin; 1388c8883902SJed Brown *B = 0; 1389d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1390c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 13917ae8954aSSatish Balay ierr = MatSetUp(A);CHKERRQ(ierr); 1392c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1393d8588912SDave May *B = A; 1394d8588912SDave May PetscFunctionReturn(0); 1395d8588912SDave May } 1396659c6bb0SJed Brown 1397629c3df2SDmitry Karpeev #undef __FUNCT__ 1398629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ" 1399cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1400629c3df2SDmitry Karpeev { 1401629c3df2SDmitry Karpeev PetscErrorCode ierr; 1402629c3df2SDmitry Karpeev Mat_Nest *nest = (Mat_Nest*)A->data; 140383b1a929SMark Adams PetscInt m,n,M,N,i,j,k,*dnnz,*onnz,rstart; 1404649b366bSFande Kong PetscInt cstart,cend; 1405629c3df2SDmitry Karpeev Mat C; 1406629c3df2SDmitry Karpeev 1407629c3df2SDmitry Karpeev PetscFunctionBegin; 1408629c3df2SDmitry Karpeev ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1409629c3df2SDmitry Karpeev ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1410649b366bSFande Kong ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); 1411629c3df2SDmitry Karpeev switch (reuse) { 1412629c3df2SDmitry Karpeev case MAT_INITIAL_MATRIX: 1413ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 1414629c3df2SDmitry Karpeev ierr = MatSetType(C,newtype);CHKERRQ(ierr); 1415629c3df2SDmitry Karpeev ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 1416629c3df2SDmitry Karpeev *newmat = C; 1417629c3df2SDmitry Karpeev break; 1418629c3df2SDmitry Karpeev case MAT_REUSE_MATRIX: 1419629c3df2SDmitry Karpeev C = *newmat; 1420629c3df2SDmitry Karpeev break; 1421ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse"); 1422629c3df2SDmitry Karpeev } 1423785e854fSJed Brown ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr); 1424629c3df2SDmitry Karpeev onnz = dnnz + m; 1425629c3df2SDmitry Karpeev for (k=0; k<m; k++) { 1426629c3df2SDmitry Karpeev dnnz[k] = 0; 1427629c3df2SDmitry Karpeev onnz[k] = 0; 1428629c3df2SDmitry Karpeev } 1429629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1430629c3df2SDmitry Karpeev IS bNis; 1431629c3df2SDmitry Karpeev PetscInt bN; 1432629c3df2SDmitry Karpeev const PetscInt *bNindices; 1433629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1434629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1435629c3df2SDmitry Karpeev ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr); 1436629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1437629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1438629c3df2SDmitry Karpeev PetscSF bmsf; 1439649b366bSFande Kong PetscSFNode *iremote; 1440629c3df2SDmitry Karpeev Mat B; 1441649b366bSFande Kong PetscInt bm, *sub_dnnz,*sub_onnz, br; 1442629c3df2SDmitry Karpeev const PetscInt *bmindices; 1443629c3df2SDmitry Karpeev B = nest->m[i][j]; 1444629c3df2SDmitry Karpeev if (!B) continue; 1445629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1446629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1447ce94432eSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr); 1448649b366bSFande Kong ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr); 1449649b366bSFande Kong ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr); 1450649b366bSFande Kong ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr); 1451649b366bSFande Kong for (k = 0; k < bm; ++k){ 1452649b366bSFande Kong sub_dnnz[k] = 0; 1453649b366bSFande Kong sub_onnz[k] = 0; 1454649b366bSFande Kong } 1455629c3df2SDmitry Karpeev /* 1456629c3df2SDmitry Karpeev Locate the owners for all of the locally-owned global row indices for this row block. 1457629c3df2SDmitry Karpeev These determine the roots of PetscSF used to communicate preallocation data to row owners. 1458629c3df2SDmitry Karpeev The roots correspond to the dnnz and onnz entries; thus, there are two roots per row. 1459629c3df2SDmitry Karpeev */ 146083b1a929SMark Adams ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 1461629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1462649b366bSFande Kong PetscInt row = bmindices[br], rowowner = 0, brncols, col; 1463629c3df2SDmitry Karpeev const PetscInt *brcols; 1464a4b3d3acSMatthew G Knepley PetscInt rowrel = 0; /* row's relative index on its owner rank */ 1465629c3df2SDmitry Karpeev ierr = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr); 1466649b366bSFande Kong /* how many roots */ 1467649b366bSFande Kong iremote[br].rank = rowowner; iremote[br].index = rowrel; /* edge from bmdnnz to dnnz */ 1468649b366bSFande Kong /* get nonzero pattern */ 146983b1a929SMark Adams ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr); 1470629c3df2SDmitry Karpeev for (k=0; k<brncols; k++) { 1471629c3df2SDmitry Karpeev col = bNindices[brcols[k]]; 1472649b366bSFande Kong if(col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]){ 1473649b366bSFande Kong sub_dnnz[br]++; 1474649b366bSFande Kong }else{ 1475649b366bSFande Kong sub_onnz[br]++; 1476649b366bSFande Kong } 1477629c3df2SDmitry Karpeev } 147883b1a929SMark Adams ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr); 1479629c3df2SDmitry Karpeev } 1480629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1481629c3df2SDmitry Karpeev /* bsf will have to take care of disposing of bedges. */ 1482649b366bSFande Kong ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1483649b366bSFande Kong ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr); 1484649b366bSFande Kong ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr); 1485649b366bSFande Kong ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr); 1486649b366bSFande Kong ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr); 1487649b366bSFande Kong ierr = PetscFree(sub_dnnz);CHKERRQ(ierr); 1488649b366bSFande Kong ierr = PetscFree(sub_onnz);CHKERRQ(ierr); 1489629c3df2SDmitry Karpeev ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr); 1490629c3df2SDmitry Karpeev } 149122d28d08SBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1492629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1493629c3df2SDmitry Karpeev } 1494629c3df2SDmitry Karpeev ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr); 1495629c3df2SDmitry Karpeev ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr); 1496629c3df2SDmitry Karpeev ierr = PetscFree(dnnz);CHKERRQ(ierr); 1497629c3df2SDmitry Karpeev 1498629c3df2SDmitry Karpeev /* Fill by row */ 1499629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1500629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1501629c3df2SDmitry Karpeev IS bNis; 1502629c3df2SDmitry Karpeev PetscInt bN; 1503629c3df2SDmitry Karpeev const PetscInt *bNindices; 1504629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1505629c3df2SDmitry Karpeev ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr); 1506629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1507629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1508629c3df2SDmitry Karpeev Mat B; 1509629c3df2SDmitry Karpeev PetscInt bm, br; 1510629c3df2SDmitry Karpeev const PetscInt *bmindices; 1511629c3df2SDmitry Karpeev B = nest->m[i][j]; 1512629c3df2SDmitry Karpeev if (!B) continue; 1513629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1514629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 151583b1a929SMark Adams ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 1516629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1517629c3df2SDmitry Karpeev PetscInt row = bmindices[br], brncols, *cols; 1518629c3df2SDmitry Karpeev const PetscInt *brcols; 1519629c3df2SDmitry Karpeev const PetscScalar *brcoldata; 152083b1a929SMark Adams ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1521785e854fSJed Brown ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr); 152226fbe8dcSKarl Rupp for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]]; 1523629c3df2SDmitry Karpeev /* 1524629c3df2SDmitry Karpeev Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match. 1525629c3df2SDmitry Karpeev Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES. 1526629c3df2SDmitry Karpeev */ 1527a2ea699eSBarry Smith ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr); 152883b1a929SMark Adams ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1529629c3df2SDmitry Karpeev ierr = PetscFree(cols);CHKERRQ(ierr); 1530629c3df2SDmitry Karpeev } 1531629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1532629c3df2SDmitry Karpeev } 1533a2ea699eSBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1534629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1535629c3df2SDmitry Karpeev } 1536629c3df2SDmitry Karpeev ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1537629c3df2SDmitry Karpeev ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1538629c3df2SDmitry Karpeev PetscFunctionReturn(0); 1539629c3df2SDmitry Karpeev } 1540629c3df2SDmitry Karpeev 1541659c6bb0SJed Brown /*MC 1542659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1543659c6bb0SJed Brown 1544659c6bb0SJed Brown Level: intermediate 1545659c6bb0SJed Brown 1546659c6bb0SJed Brown Notes: 1547659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1548659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1549950540a4SJed Brown It is usually used with DMComposite and DMCreateMatrix() 1550659c6bb0SJed Brown 1551659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1552659c6bb0SJed Brown M*/ 1553c8883902SJed Brown #undef __FUNCT__ 1554c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 15558cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A) 1556c8883902SJed Brown { 1557c8883902SJed Brown Mat_Nest *s; 1558c8883902SJed Brown PetscErrorCode ierr; 1559c8883902SJed Brown 1560c8883902SJed Brown PetscFunctionBegin; 1561b00a9115SJed Brown ierr = PetscNewLog(A,&s);CHKERRQ(ierr); 1562c8883902SJed Brown A->data = (void*)s; 1563e7c19651SJed Brown 1564e7c19651SJed Brown s->nr = -1; 1565e7c19651SJed Brown s->nc = -1; 15660298fd71SBarry Smith s->m = NULL; 1567e7c19651SJed Brown s->splitassembly = PETSC_FALSE; 1568c8883902SJed Brown 1569c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 157026fbe8dcSKarl Rupp 1571c8883902SJed Brown A->ops->mult = MatMult_Nest; 15729194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 1573c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 15749194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1575c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1576c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1577c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1578c222c20dSDavid Ham A->ops->copy = MatCopy_Nest; 1579c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1580c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1581c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1582c8883902SJed Brown A->ops->view = MatView_Nest; 1583c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1584c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1585c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1586429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 1587429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 1588a061e289SJed Brown A->ops->scale = MatScale_Nest; 1589a061e289SJed Brown A->ops->shift = MatShift_Nest; 1590c8883902SJed Brown 1591c8883902SJed Brown A->spptr = 0; 1592c8883902SJed Brown A->assembled = PETSC_FALSE; 1593c8883902SJed Brown 1594c8883902SJed Brown /* expose Nest api's */ 1595bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C", MatNestGetSubMat_Nest);CHKERRQ(ierr); 1596bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1597bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1598bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C", MatNestGetSize_Nest);CHKERRQ(ierr); 1599bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C", MatNestGetISs_Nest);CHKERRQ(ierr); 1600bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 1601bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C", MatNestSetVecType_Nest);CHKERRQ(ierr); 1602bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C", MatNestSetSubMats_Nest);CHKERRQ(ierr); 160383b1a929SMark Adams ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr); 1604c8883902SJed Brown 1605c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1606c8883902SJed Brown PetscFunctionReturn(0); 1607c8883902SJed Brown } 1608