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[]); 67874fa86SDave May static PetscErrorCode MatGetVecs_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 } 364ce94432eSBarry Smith ierr = MPI_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr); 3658188e55aSJed Brown 366e072481dSJed Brown if (isFullColGlobal) { 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); 373f349c1fdSJed Brown *B = vs->m[row][col]; 3748188e55aSJed Brown } 375f349c1fdSJed Brown PetscFunctionReturn(0); 376f349c1fdSJed Brown } 377f349c1fdSJed Brown 378f349c1fdSJed Brown #undef __FUNCT__ 379f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 380f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 381f349c1fdSJed Brown { 382f349c1fdSJed Brown PetscErrorCode ierr; 383f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 384f349c1fdSJed Brown Mat sub; 385f349c1fdSJed Brown 386f349c1fdSJed Brown PetscFunctionBegin; 387f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 388f349c1fdSJed Brown switch (reuse) { 389f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 3907874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 391f349c1fdSJed Brown *B = sub; 392f349c1fdSJed Brown break; 393f349c1fdSJed Brown case MAT_REUSE_MATRIX: 394ce94432eSBarry Smith if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 395f349c1fdSJed Brown break; 396f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 397f349c1fdSJed Brown break; 398f349c1fdSJed Brown } 399f349c1fdSJed Brown PetscFunctionReturn(0); 400f349c1fdSJed Brown } 401f349c1fdSJed Brown 402f349c1fdSJed Brown #undef __FUNCT__ 403f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 404f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 405f349c1fdSJed Brown { 406f349c1fdSJed Brown PetscErrorCode ierr; 407f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 408f349c1fdSJed Brown Mat sub; 409f349c1fdSJed Brown 410f349c1fdSJed Brown PetscFunctionBegin; 411f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 412f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 413f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 414f349c1fdSJed Brown *B = sub; 415d8588912SDave May PetscFunctionReturn(0); 416d8588912SDave May } 417d8588912SDave May 418d8588912SDave May #undef __FUNCT__ 419d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 420207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 421d8588912SDave May { 422d8588912SDave May PetscErrorCode ierr; 423f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 424f349c1fdSJed Brown Mat sub; 425d8588912SDave May 426d8588912SDave May PetscFunctionBegin; 427f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 428ce94432eSBarry Smith if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 429f349c1fdSJed Brown if (sub) { 430ce94432eSBarry Smith if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4316bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 432d8588912SDave May } 433d8588912SDave May PetscFunctionReturn(0); 434d8588912SDave May } 435d8588912SDave May 436d8588912SDave May #undef __FUNCT__ 4377874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 4387874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 4397874fa86SDave May { 4407874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4417874fa86SDave May PetscInt i; 4427874fa86SDave May PetscErrorCode ierr; 4437874fa86SDave May 4447874fa86SDave May PetscFunctionBegin; 4457874fa86SDave May for (i=0; i<bA->nr; i++) { 446429bac76SJed Brown Vec bv; 447429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4487874fa86SDave May if (bA->m[i][i]) { 449429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 4507874fa86SDave May } else { 451429bac76SJed Brown ierr = VecSet(bv,1.0);CHKERRQ(ierr); 4527874fa86SDave May } 453429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4547874fa86SDave May } 4557874fa86SDave May PetscFunctionReturn(0); 4567874fa86SDave May } 4577874fa86SDave May 4587874fa86SDave May #undef __FUNCT__ 4597874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 4607874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 4617874fa86SDave May { 4627874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 463429bac76SJed Brown Vec bl,*br; 4647874fa86SDave May PetscInt i,j; 4657874fa86SDave May PetscErrorCode ierr; 4667874fa86SDave May 4677874fa86SDave May PetscFunctionBegin; 4683f800ebeSJed Brown ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr); 4692e6472ebSElliott Sales de Andrade if (r) { 470429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 4712e6472ebSElliott Sales de Andrade } 4722e6472ebSElliott Sales de Andrade bl = NULL; 4737874fa86SDave May for (i=0; i<bA->nr; i++) { 4742e6472ebSElliott Sales de Andrade if (l) { 475429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4762e6472ebSElliott Sales de Andrade } 4777874fa86SDave May for (j=0; j<bA->nc; j++) { 4787874fa86SDave May if (bA->m[i][j]) { 479429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 4807874fa86SDave May } 4817874fa86SDave May } 4822e6472ebSElliott Sales de Andrade if (l) { 483a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4847874fa86SDave May } 4852e6472ebSElliott Sales de Andrade } 4862e6472ebSElliott Sales de Andrade if (r) { 487429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 4882e6472ebSElliott Sales de Andrade } 489429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 4907874fa86SDave May PetscFunctionReturn(0); 4917874fa86SDave May } 4927874fa86SDave May 4937874fa86SDave May #undef __FUNCT__ 494a061e289SJed Brown #define __FUNCT__ "MatScale_Nest" 495a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 496a061e289SJed Brown { 497a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 498a061e289SJed Brown PetscInt i,j; 499a061e289SJed Brown PetscErrorCode ierr; 500a061e289SJed Brown 501a061e289SJed Brown PetscFunctionBegin; 502a061e289SJed Brown for (i=0; i<bA->nr; i++) { 503a061e289SJed Brown for (j=0; j<bA->nc; j++) { 504a061e289SJed Brown if (bA->m[i][j]) { 505a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 506a061e289SJed Brown } 507a061e289SJed Brown } 508a061e289SJed Brown } 509a061e289SJed Brown PetscFunctionReturn(0); 510a061e289SJed Brown } 511a061e289SJed Brown 512a061e289SJed Brown #undef __FUNCT__ 513a061e289SJed Brown #define __FUNCT__ "MatShift_Nest" 514a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 515a061e289SJed Brown { 516a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 517a061e289SJed Brown PetscInt i; 518a061e289SJed Brown PetscErrorCode ierr; 519a061e289SJed Brown 520a061e289SJed Brown PetscFunctionBegin; 521a061e289SJed Brown for (i=0; i<bA->nr; i++) { 522ce94432eSBarry 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); 523a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 524a061e289SJed Brown } 525a061e289SJed Brown PetscFunctionReturn(0); 526a061e289SJed Brown } 527a061e289SJed Brown 528a061e289SJed Brown #undef __FUNCT__ 529d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest" 530207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 531d8588912SDave May { 532d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 533d8588912SDave May Vec *L,*R; 534d8588912SDave May MPI_Comm comm; 535d8588912SDave May PetscInt i,j; 536d8588912SDave May PetscErrorCode ierr; 537d8588912SDave May 538d8588912SDave May PetscFunctionBegin; 539ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 540d8588912SDave May if (right) { 541d8588912SDave May /* allocate R */ 542d8588912SDave May ierr = PetscMalloc(sizeof(Vec) * bA->nc, &R);CHKERRQ(ierr); 543d8588912SDave May /* Create the right vectors */ 544d8588912SDave May for (j=0; j<bA->nc; j++) { 545d8588912SDave May for (i=0; i<bA->nr; i++) { 546d8588912SDave May if (bA->m[i][j]) { 5470298fd71SBarry Smith ierr = MatGetVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr); 548d8588912SDave May break; 549d8588912SDave May } 550d8588912SDave May } 551d8588912SDave May if (i==bA->nr) { 552d8588912SDave May /* have an empty column */ 553ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 554d8588912SDave May } 555d8588912SDave May } 556f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 557d8588912SDave May /* hand back control to the nest vector */ 558d8588912SDave May for (j=0; j<bA->nc; j++) { 5596bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 560d8588912SDave May } 561d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 562d8588912SDave May } 563d8588912SDave May 564d8588912SDave May if (left) { 565d8588912SDave May /* allocate L */ 566d8588912SDave May ierr = PetscMalloc(sizeof(Vec) * bA->nr, &L);CHKERRQ(ierr); 567d8588912SDave May /* Create the left vectors */ 568d8588912SDave May for (i=0; i<bA->nr; i++) { 569d8588912SDave May for (j=0; j<bA->nc; j++) { 570d8588912SDave May if (bA->m[i][j]) { 5710298fd71SBarry Smith ierr = MatGetVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr); 572d8588912SDave May break; 573d8588912SDave May } 574d8588912SDave May } 575d8588912SDave May if (j==bA->nc) { 576d8588912SDave May /* have an empty row */ 577ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 578d8588912SDave May } 579d8588912SDave May } 580d8588912SDave May 581f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 582d8588912SDave May for (i=0; i<bA->nr; i++) { 5836bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 584d8588912SDave May } 585d8588912SDave May 586d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 587d8588912SDave May } 588d8588912SDave May PetscFunctionReturn(0); 589d8588912SDave May } 590d8588912SDave May 591d8588912SDave May #undef __FUNCT__ 592d8588912SDave May #define __FUNCT__ "MatView_Nest" 593207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 594d8588912SDave May { 595d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 596d8588912SDave May PetscBool isascii; 597d8588912SDave May PetscInt i,j; 598d8588912SDave May PetscErrorCode ierr; 599d8588912SDave May 600d8588912SDave May PetscFunctionBegin; 601251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 602d8588912SDave May if (isascii) { 603d8588912SDave May 604d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n"); 605d8588912SDave May PetscViewerASCIIPushTab(viewer); /* push0 */ 606d8588912SDave May PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 607d8588912SDave May 608d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n"); 609d8588912SDave May for (i=0; i<bA->nr; i++) { 610d8588912SDave May for (j=0; j<bA->nc; j++) { 61119fd82e9SBarry Smith MatType type; 612270f95d7SJed Brown char name[256] = "",prefix[256] = ""; 613d8588912SDave May PetscInt NR,NC; 614d8588912SDave May PetscBool isNest = PETSC_FALSE; 615d8588912SDave May 616d8588912SDave May if (!bA->m[i][j]) { 6170298fd71SBarry Smith PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j); 618d8588912SDave May continue; 619d8588912SDave May } 620d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 621d8588912SDave May ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr); 6228caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 6238caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 624251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 625d8588912SDave May 626270f95d7SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 627d8588912SDave May 628d8588912SDave May if (isNest) { 629270f95d7SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 630d8588912SDave May ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 631270f95d7SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 632d8588912SDave May } 633d8588912SDave May } 634d8588912SDave May } 635d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 636d8588912SDave May } 637d8588912SDave May PetscFunctionReturn(0); 638d8588912SDave May } 639d8588912SDave May 640d8588912SDave May #undef __FUNCT__ 641d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 642207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 643d8588912SDave May { 644d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 645d8588912SDave May PetscInt i,j; 646d8588912SDave May PetscErrorCode ierr; 647d8588912SDave May 648d8588912SDave May PetscFunctionBegin; 649d8588912SDave May for (i=0; i<bA->nr; i++) { 650d8588912SDave May for (j=0; j<bA->nc; j++) { 651d8588912SDave May if (!bA->m[i][j]) continue; 652d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 653d8588912SDave May } 654d8588912SDave May } 655d8588912SDave May PetscFunctionReturn(0); 656d8588912SDave May } 657d8588912SDave May 658d8588912SDave May #undef __FUNCT__ 659*c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest" 660*c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str) 661*c222c20dSDavid Ham { 662*c222c20dSDavid Ham Mat_Nest *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data; 663*c222c20dSDavid Ham PetscInt i,j,nr = bA->nr,nc = bA->nc; 664*c222c20dSDavid Ham PetscErrorCode ierr; 665*c222c20dSDavid Ham 666*c222c20dSDavid Ham PetscFunctionBegin; 667*c222c20dSDavid 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); 668*c222c20dSDavid Ham for (i=0; i<nr; i++) { 669*c222c20dSDavid Ham for (j=0; j<nc; j++) { 670*c222c20dSDavid Ham ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr); 671*c222c20dSDavid Ham } 672*c222c20dSDavid Ham } 673*c222c20dSDavid Ham PetscFunctionReturn(0); 674*c222c20dSDavid Ham } 675*c222c20dSDavid Ham 676*c222c20dSDavid Ham #undef __FUNCT__ 677d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 678207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 679d8588912SDave May { 680d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 681841e96a3SJed Brown Mat *b; 682841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 683d8588912SDave May PetscErrorCode ierr; 684d8588912SDave May 685d8588912SDave May PetscFunctionBegin; 686785e854fSJed Brown ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr); 687841e96a3SJed Brown for (i=0; i<nr; i++) { 688841e96a3SJed Brown for (j=0; j<nc; j++) { 689841e96a3SJed Brown if (bA->m[i][j]) { 690841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 691841e96a3SJed Brown } else { 6920298fd71SBarry Smith b[i*nc+j] = NULL; 693d8588912SDave May } 694d8588912SDave May } 695d8588912SDave May } 696ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 697841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 698841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 6996bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 700d8588912SDave May } 701d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 702d8588912SDave May 703841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 704841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 705d8588912SDave May PetscFunctionReturn(0); 706d8588912SDave May } 707d8588912SDave May 708d8588912SDave May /* nest api */ 709d8588912SDave May #undef __FUNCT__ 710d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 711d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 712d8588912SDave May { 713d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 7145fd66863SKarl Rupp 715d8588912SDave May PetscFunctionBegin; 716ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 717ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 718d8588912SDave May *mat = bA->m[idxm][jdxm]; 719d8588912SDave May PetscFunctionReturn(0); 720d8588912SDave May } 721d8588912SDave May 722d8588912SDave May #undef __FUNCT__ 723d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 7249ba0d327SJed Brown /*@ 725d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 726d8588912SDave May 727d8588912SDave May Not collective 728d8588912SDave May 729d8588912SDave May Input Parameters: 730629881c0SJed Brown + A - nest matrix 731d8588912SDave May . idxm - index of the matrix within the nest matrix 732629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 733d8588912SDave May 734d8588912SDave May Output Parameter: 735d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 736d8588912SDave May 737d8588912SDave May Level: developer 738d8588912SDave May 739d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 740d8588912SDave May @*/ 7417087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 742d8588912SDave May { 743699a902aSJed Brown PetscErrorCode ierr; 744d8588912SDave May 745d8588912SDave May PetscFunctionBegin; 746699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 747d8588912SDave May PetscFunctionReturn(0); 748d8588912SDave May } 749d8588912SDave May 750d8588912SDave May #undef __FUNCT__ 7510782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 7520782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 7530782ca92SJed Brown { 7540782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 7550782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 7560782ca92SJed Brown PetscErrorCode ierr; 7570782ca92SJed Brown 7580782ca92SJed Brown PetscFunctionBegin; 759ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 760ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 7610782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 7620782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 7630782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 7640782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 7650782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 7660782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 767ce94432eSBarry 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); 768ce94432eSBarry 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); 76926fbe8dcSKarl Rupp 7700782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 7710782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 7720782ca92SJed Brown bA->m[idxm][jdxm] = mat; 7730782ca92SJed Brown PetscFunctionReturn(0); 7740782ca92SJed Brown } 7750782ca92SJed Brown 7760782ca92SJed Brown #undef __FUNCT__ 7770782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 7789ba0d327SJed Brown /*@ 7790782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 7800782ca92SJed Brown 7810782ca92SJed Brown Logically collective on the submatrix communicator 7820782ca92SJed Brown 7830782ca92SJed Brown Input Parameters: 7840782ca92SJed Brown + A - nest matrix 7850782ca92SJed Brown . idxm - index of the matrix within the nest matrix 7860782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 7870782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 7880782ca92SJed Brown 7890782ca92SJed Brown Notes: 7900782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 7910782ca92SJed Brown 7920782ca92SJed Brown This increments the reference count of the submatrix. 7930782ca92SJed Brown 7940782ca92SJed Brown Level: developer 7950782ca92SJed Brown 7960782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 7970782ca92SJed Brown @*/ 7980782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 7990782ca92SJed Brown { 8000782ca92SJed Brown PetscErrorCode ierr; 8010782ca92SJed Brown 8020782ca92SJed Brown PetscFunctionBegin; 8030782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 8040782ca92SJed Brown PetscFunctionReturn(0); 8050782ca92SJed Brown } 8060782ca92SJed Brown 8070782ca92SJed Brown #undef __FUNCT__ 808d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 809d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 810d8588912SDave May { 811d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 8125fd66863SKarl Rupp 813d8588912SDave May PetscFunctionBegin; 81426fbe8dcSKarl Rupp if (M) *M = bA->nr; 81526fbe8dcSKarl Rupp if (N) *N = bA->nc; 81626fbe8dcSKarl Rupp if (mat) *mat = bA->m; 817d8588912SDave May PetscFunctionReturn(0); 818d8588912SDave May } 819d8588912SDave May 820d8588912SDave May #undef __FUNCT__ 821d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 822d8588912SDave May /*@C 823d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 824d8588912SDave May 825d8588912SDave May Not collective 826d8588912SDave May 827d8588912SDave May Input Parameters: 828629881c0SJed Brown . A - nest matrix 829d8588912SDave May 830d8588912SDave May Output Parameter: 831629881c0SJed Brown + M - number of rows in the nest matrix 832d8588912SDave May . N - number of cols in the nest matrix 833629881c0SJed Brown - mat - 2d array of matrices 834d8588912SDave May 835d8588912SDave May Notes: 836d8588912SDave May 837d8588912SDave May The user should not free the array mat. 838d8588912SDave May 839d8588912SDave May Level: developer 840d8588912SDave May 841d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 842d8588912SDave May @*/ 8437087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 844d8588912SDave May { 845699a902aSJed Brown PetscErrorCode ierr; 846d8588912SDave May 847d8588912SDave May PetscFunctionBegin; 848699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 849d8588912SDave May PetscFunctionReturn(0); 850d8588912SDave May } 851d8588912SDave May 852d8588912SDave May #undef __FUNCT__ 853d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 8547087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 855d8588912SDave May { 856d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 857d8588912SDave May 858d8588912SDave May PetscFunctionBegin; 85926fbe8dcSKarl Rupp if (M) *M = bA->nr; 86026fbe8dcSKarl Rupp if (N) *N = bA->nc; 861d8588912SDave May PetscFunctionReturn(0); 862d8588912SDave May } 863d8588912SDave May 864d8588912SDave May #undef __FUNCT__ 865d8588912SDave May #define __FUNCT__ "MatNestGetSize" 8669ba0d327SJed Brown /*@ 867d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 868d8588912SDave May 869d8588912SDave May Not collective 870d8588912SDave May 871d8588912SDave May Input Parameters: 872d8588912SDave May . A - nest matrix 873d8588912SDave May 874d8588912SDave May Output Parameter: 875629881c0SJed Brown + M - number of rows in the nested mat 876629881c0SJed Brown - N - number of cols in the nested mat 877d8588912SDave May 878d8588912SDave May Notes: 879d8588912SDave May 880d8588912SDave May Level: developer 881d8588912SDave May 882d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 883d8588912SDave May @*/ 8847087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 885d8588912SDave May { 886699a902aSJed Brown PetscErrorCode ierr; 887d8588912SDave May 888d8588912SDave May PetscFunctionBegin; 889699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 890d8588912SDave May PetscFunctionReturn(0); 891d8588912SDave May } 892d8588912SDave May 893900e7ff2SJed Brown #undef __FUNCT__ 894900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest" 895f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 896900e7ff2SJed Brown { 897900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 898900e7ff2SJed Brown PetscInt i; 899900e7ff2SJed Brown 900900e7ff2SJed Brown PetscFunctionBegin; 901900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 902900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 903900e7ff2SJed Brown PetscFunctionReturn(0); 904900e7ff2SJed Brown } 905900e7ff2SJed Brown 906900e7ff2SJed Brown #undef __FUNCT__ 907900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs" 908900e7ff2SJed Brown /*@C 909900e7ff2SJed Brown MatNestGetISs - Returns the index sets partitioning the row and column spaces 910900e7ff2SJed Brown 911900e7ff2SJed Brown Not collective 912900e7ff2SJed Brown 913900e7ff2SJed Brown Input Parameters: 914900e7ff2SJed Brown . A - nest matrix 915900e7ff2SJed Brown 916900e7ff2SJed Brown Output Parameter: 917900e7ff2SJed Brown + rows - array of row index sets 918900e7ff2SJed Brown - cols - array of column index sets 919900e7ff2SJed Brown 920900e7ff2SJed Brown Level: advanced 921900e7ff2SJed Brown 922900e7ff2SJed Brown Notes: 923900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 924900e7ff2SJed Brown 925900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs() 926900e7ff2SJed Brown @*/ 927900e7ff2SJed Brown PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 928900e7ff2SJed Brown { 929900e7ff2SJed Brown PetscErrorCode ierr; 930900e7ff2SJed Brown 931900e7ff2SJed Brown PetscFunctionBegin; 932900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 933900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 934900e7ff2SJed Brown PetscFunctionReturn(0); 935900e7ff2SJed Brown } 936900e7ff2SJed Brown 937900e7ff2SJed Brown #undef __FUNCT__ 938900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest" 939f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 940900e7ff2SJed Brown { 941900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 942900e7ff2SJed Brown PetscInt i; 943900e7ff2SJed Brown 944900e7ff2SJed Brown PetscFunctionBegin; 945900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 946900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 947900e7ff2SJed Brown PetscFunctionReturn(0); 948900e7ff2SJed Brown } 949900e7ff2SJed Brown 950900e7ff2SJed Brown #undef __FUNCT__ 951900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs" 952900e7ff2SJed Brown /*@C 953900e7ff2SJed Brown MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 954900e7ff2SJed Brown 955900e7ff2SJed Brown Not collective 956900e7ff2SJed Brown 957900e7ff2SJed Brown Input Parameters: 958900e7ff2SJed Brown . A - nest matrix 959900e7ff2SJed Brown 960900e7ff2SJed Brown Output Parameter: 9610298fd71SBarry Smith + rows - array of row index sets (or NULL to ignore) 9620298fd71SBarry Smith - cols - array of column index sets (or NULL to ignore) 963900e7ff2SJed Brown 964900e7ff2SJed Brown Level: advanced 965900e7ff2SJed Brown 966900e7ff2SJed Brown Notes: 967900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 968900e7ff2SJed Brown 969900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs() 970900e7ff2SJed Brown @*/ 971900e7ff2SJed Brown PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 972900e7ff2SJed Brown { 973900e7ff2SJed Brown PetscErrorCode ierr; 974900e7ff2SJed Brown 975900e7ff2SJed Brown PetscFunctionBegin; 976900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 977900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 978900e7ff2SJed Brown PetscFunctionReturn(0); 979900e7ff2SJed Brown } 980900e7ff2SJed Brown 981207556f9SJed Brown #undef __FUNCT__ 982207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 98319fd82e9SBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,VecType vtype) 984207556f9SJed Brown { 985207556f9SJed Brown PetscErrorCode ierr; 986207556f9SJed Brown PetscBool flg; 987207556f9SJed Brown 988207556f9SJed Brown PetscFunctionBegin; 989207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 990207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 99112b53f24SSatish Balay if (flg) A->ops->getvecs = MatGetVecs_Nest; 99212b53f24SSatish Balay else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0; 993207556f9SJed Brown PetscFunctionReturn(0); 994207556f9SJed Brown } 995207556f9SJed Brown 996207556f9SJed Brown #undef __FUNCT__ 997207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 998207556f9SJed Brown /*@C 999207556f9SJed Brown MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 1000207556f9SJed Brown 1001207556f9SJed Brown Not collective 1002207556f9SJed Brown 1003207556f9SJed Brown Input Parameters: 1004207556f9SJed Brown + A - nest matrix 1005207556f9SJed Brown - vtype - type to use for creating vectors 1006207556f9SJed Brown 1007207556f9SJed Brown Notes: 1008207556f9SJed Brown 1009207556f9SJed Brown Level: developer 1010207556f9SJed Brown 1011207556f9SJed Brown .seealso: MatGetVecs() 1012207556f9SJed Brown @*/ 101319fd82e9SBarry Smith PetscErrorCode MatNestSetVecType(Mat A,VecType vtype) 1014207556f9SJed Brown { 1015207556f9SJed Brown PetscErrorCode ierr; 1016207556f9SJed Brown 1017207556f9SJed Brown PetscFunctionBegin; 101819fd82e9SBarry Smith ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr); 1019207556f9SJed Brown PetscFunctionReturn(0); 1020207556f9SJed Brown } 1021207556f9SJed Brown 1022d8588912SDave May #undef __FUNCT__ 1023c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 1024c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1025d8588912SDave May { 1026c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 1027c8883902SJed Brown PetscInt i,j,m,n,M,N; 1028d8588912SDave May PetscErrorCode ierr; 1029d8588912SDave May 1030d8588912SDave May PetscFunctionBegin; 1031c8883902SJed Brown s->nr = nr; 1032c8883902SJed Brown s->nc = nc; 1033d8588912SDave May 1034c8883902SJed Brown /* Create space for submatrices */ 1035c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 1036c8883902SJed Brown for (i=0; i<nr; i++) { 1037c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 1038d8588912SDave May } 1039c8883902SJed Brown for (i=0; i<nr; i++) { 1040c8883902SJed Brown for (j=0; j<nc; j++) { 1041c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1042c8883902SJed Brown if (a[i*nc+j]) { 1043c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1044d8588912SDave May } 1045d8588912SDave May } 1046d8588912SDave May } 1047d8588912SDave May 10488188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1049d8588912SDave May 1050c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 1051c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 1052c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1053c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1054d8588912SDave May 10558188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1056d8588912SDave May 1057c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1058c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1059c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1060c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1061c8883902SJed Brown 1062c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1063c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1064c8883902SJed Brown 10651795a4d1SJed Brown ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr); 1066d8588912SDave May PetscFunctionReturn(0); 1067d8588912SDave May } 1068d8588912SDave May 1069c8883902SJed Brown #undef __FUNCT__ 1070c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 1071c8883902SJed Brown /*@ 1072c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1073c8883902SJed Brown 1074c8883902SJed Brown Collective on Mat 1075c8883902SJed Brown 1076c8883902SJed Brown Input Parameter: 1077c8883902SJed Brown + N - nested matrix 1078c8883902SJed Brown . nr - number of nested row blocks 10790298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1080c8883902SJed Brown . nc - number of nested column blocks 10810298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 10820298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1083c8883902SJed Brown 1084c8883902SJed Brown Level: advanced 1085c8883902SJed Brown 1086c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 1087c8883902SJed Brown @*/ 1088c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1089c8883902SJed Brown { 1090c8883902SJed Brown PetscErrorCode ierr; 1091c8883902SJed Brown PetscInt i; 1092c8883902SJed Brown 1093c8883902SJed Brown PetscFunctionBegin; 1094c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1095ce94432eSBarry Smith if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1096c8883902SJed Brown if (nr && is_row) { 1097c8883902SJed Brown PetscValidPointer(is_row,3); 1098c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1099c8883902SJed Brown } 1100ce94432eSBarry Smith if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 11011664e352SJed Brown if (nc && is_col) { 1102c8883902SJed Brown PetscValidPointer(is_col,5); 1103c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1104c8883902SJed Brown } 1105c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 1106c8883902SJed 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); 1107c8883902SJed Brown PetscFunctionReturn(0); 1108c8883902SJed Brown } 1109d8588912SDave May 111077019fcaSJed Brown #undef __FUNCT__ 111177019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private" 111277019fcaSJed Brown static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb) 111377019fcaSJed Brown { 111477019fcaSJed Brown PetscErrorCode ierr; 111577019fcaSJed Brown PetscBool flg; 111677019fcaSJed Brown PetscInt i,j,m,mi,*ix; 111777019fcaSJed Brown 111877019fcaSJed Brown PetscFunctionBegin; 111977019fcaSJed Brown for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 112077019fcaSJed Brown if (islocal[i]) { 112177019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 112277019fcaSJed Brown flg = PETSC_TRUE; /* We found a non-trivial entry */ 112377019fcaSJed Brown } else { 112477019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 112577019fcaSJed Brown } 112677019fcaSJed Brown m += mi; 112777019fcaSJed Brown } 112877019fcaSJed Brown if (flg) { 1129785e854fSJed Brown ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr); 113077019fcaSJed Brown for (i=0,n=0; i<n; i++) { 11310298fd71SBarry Smith ISLocalToGlobalMapping smap = NULL; 113277019fcaSJed Brown VecScatter scat; 113377019fcaSJed Brown IS isreq; 113477019fcaSJed Brown Vec lvec,gvec; 11353361c9a7SJed Brown union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x; 113677019fcaSJed Brown Mat sub; 113777019fcaSJed Brown 1138ce94432eSBarry Smith if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers"); 113977019fcaSJed Brown if (colflg) { 114077019fcaSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 114177019fcaSJed Brown } else { 114277019fcaSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 114377019fcaSJed Brown } 11440298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);} 114577019fcaSJed Brown if (islocal[i]) { 114677019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 114777019fcaSJed Brown } else { 114877019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 114977019fcaSJed Brown } 115077019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = j; 115177019fcaSJed Brown if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);} 115277019fcaSJed Brown /* 115377019fcaSJed Brown Now we need to extract the monolithic global indices that correspond to the given split global indices. 115477019fcaSJed 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. 115577019fcaSJed Brown The approach here is ugly because it uses VecScatter to move indices. 115677019fcaSJed Brown */ 115777019fcaSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr); 115877019fcaSJed Brown ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr); 115977019fcaSJed Brown ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr); 11600298fd71SBarry Smith ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr); 116177019fcaSJed Brown ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 116277019fcaSJed Brown for (j=0; j<mi; j++) x[j].integer = ix[m+j]; 116377019fcaSJed Brown ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 116477019fcaSJed Brown ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 116577019fcaSJed Brown ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 116677019fcaSJed Brown ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 116777019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = x[j].integer; 116877019fcaSJed Brown ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 116977019fcaSJed Brown ierr = VecDestroy(&lvec);CHKERRQ(ierr); 117077019fcaSJed Brown ierr = VecDestroy(&gvec);CHKERRQ(ierr); 117177019fcaSJed Brown ierr = ISDestroy(&isreq);CHKERRQ(ierr); 117277019fcaSJed Brown ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 117377019fcaSJed Brown m += mi; 117477019fcaSJed Brown } 1175ce94432eSBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 11760298fd71SBarry Smith *ltogb = NULL; 117777019fcaSJed Brown } else { 11780298fd71SBarry Smith *ltog = NULL; 11790298fd71SBarry Smith *ltogb = NULL; 118077019fcaSJed Brown } 118177019fcaSJed Brown PetscFunctionReturn(0); 118277019fcaSJed Brown } 118377019fcaSJed Brown 118477019fcaSJed Brown 1185d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1186d8588912SDave May /* 1187d8588912SDave May nprocessors = NP 1188d8588912SDave May Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1)) 1189d8588912SDave May proc 0: => (g_0,h_0,) 1190d8588912SDave May proc 1: => (g_1,h_1,) 1191d8588912SDave May ... 1192d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1193d8588912SDave May 1194d8588912SDave May proc 0: proc 1: proc nprocs-1: 1195d8588912SDave May is[0] = (0,1,2,...,nlocal(g_0)-1) (0,1,...,nlocal(g_1)-1) (0,1,...,nlocal(g_NP-1)) 1196d8588912SDave May 1197d8588912SDave May proc 0: 1198d8588912SDave May is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1) 1199d8588912SDave May proc 1: 1200d8588912SDave May is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1) 1201d8588912SDave May 1202d8588912SDave May proc NP-1: 1203d8588912SDave May is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1) 1204d8588912SDave May */ 1205d8588912SDave May #undef __FUNCT__ 1206d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1207841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1208d8588912SDave May { 1209e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 12108188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1211d8588912SDave May PetscErrorCode ierr; 12120298fd71SBarry Smith Mat sub = NULL; 1213d8588912SDave May 1214d8588912SDave May PetscFunctionBegin; 12158188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 12168188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1217d8588912SDave May if (is_row) { /* valid IS is passed in */ 1218d8588912SDave May /* refs on is[] are incremeneted */ 1219e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1220d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 122126fbe8dcSKarl Rupp 1222e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1223d8588912SDave May } 12242ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 12258188e55aSJed Brown nsum = 0; 12268188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 12278188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1228ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 12290298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 1230ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12318188e55aSJed Brown nsum += n; 12328188e55aSJed Brown } 1233ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 123430bc264bSJed Brown offset -= nsum; 1235e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1236f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12370298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 12382ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1239ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1240e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 12412ae74bdbSJed Brown offset += n; 1242d8588912SDave May } 1243d8588912SDave May } 1244d8588912SDave May 1245d8588912SDave May if (is_col) { /* valid IS is passed in */ 1246d8588912SDave May /* refs on is[] are incremeneted */ 1247e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1248d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 124926fbe8dcSKarl Rupp 1250e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1251d8588912SDave May } 12522ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 12532ae74bdbSJed Brown offset = A->cmap->rstart; 12548188e55aSJed Brown nsum = 0; 12558188e55aSJed Brown for (j=0; j<vs->nc; j++) { 12568188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1257ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 12580298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 1259ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12608188e55aSJed Brown nsum += n; 12618188e55aSJed Brown } 1262ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 126330bc264bSJed Brown offset -= nsum; 1264e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1265f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 12660298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 12672ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1268ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1269e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 12702ae74bdbSJed Brown offset += n; 1271d8588912SDave May } 1272d8588912SDave May } 1273e2d7f03fSJed Brown 1274e2d7f03fSJed Brown /* Set up the local ISs */ 1275785e854fSJed Brown ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 1276785e854fSJed Brown ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 1277e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1278e2d7f03fSJed Brown IS isloc; 12790298fd71SBarry Smith ISLocalToGlobalMapping rmap = NULL; 1280e2d7f03fSJed Brown PetscInt nlocal,bs; 1281e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12820298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);} 1283207556f9SJed Brown if (rmap) { 1284e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1285e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1286e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1287e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1288207556f9SJed Brown } else { 1289207556f9SJed Brown nlocal = 0; 12900298fd71SBarry Smith isloc = NULL; 1291207556f9SJed Brown } 1292e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1293e2d7f03fSJed Brown offset += nlocal; 1294e2d7f03fSJed Brown } 12958188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1296e2d7f03fSJed Brown IS isloc; 12970298fd71SBarry Smith ISLocalToGlobalMapping cmap = NULL; 1298e2d7f03fSJed Brown PetscInt nlocal,bs; 1299e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 13000298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);} 1301207556f9SJed Brown if (cmap) { 1302e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1303e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1304e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1305e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1306207556f9SJed Brown } else { 1307207556f9SJed Brown nlocal = 0; 13080298fd71SBarry Smith isloc = NULL; 1309207556f9SJed Brown } 1310e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1311e2d7f03fSJed Brown offset += nlocal; 1312e2d7f03fSJed Brown } 13130189643fSJed Brown 131477019fcaSJed Brown /* Set up the aggregate ISLocalToGlobalMapping */ 131577019fcaSJed Brown { 131677019fcaSJed Brown ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb; 131777019fcaSJed Brown ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr); 131877019fcaSJed Brown ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr); 131977019fcaSJed Brown if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 132077019fcaSJed Brown if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);} 132177019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 132277019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr); 132377019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 132477019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr); 132577019fcaSJed Brown } 132677019fcaSJed Brown 13270189643fSJed Brown #if defined(PETSC_USE_DEBUG) 13280189643fSJed Brown for (i=0; i<vs->nr; i++) { 13290189643fSJed Brown for (j=0; j<vs->nc; j++) { 13300189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 13310189643fSJed Brown Mat B = vs->m[i][j]; 13320189643fSJed Brown if (!B) continue; 13330189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 13340189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 13350189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 13360189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 13370189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 13380189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 1339ce94432eSBarry 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); 1340ce94432eSBarry 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); 13410189643fSJed Brown } 13420189643fSJed Brown } 13430189643fSJed Brown #endif 1344a061e289SJed Brown 1345a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1346a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1347a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1348a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1349a061e289SJed Brown } 1350a061e289SJed Brown } 1351a061e289SJed Brown A->assembled = PETSC_TRUE; 1352d8588912SDave May PetscFunctionReturn(0); 1353d8588912SDave May } 1354d8588912SDave May 1355d8588912SDave May #undef __FUNCT__ 1356d8588912SDave May #define __FUNCT__ "MatCreateNest" 135745c38901SJed Brown /*@C 1358659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1359659c6bb0SJed Brown 1360659c6bb0SJed Brown Collective on Mat 1361659c6bb0SJed Brown 1362659c6bb0SJed Brown Input Parameter: 1363659c6bb0SJed Brown + comm - Communicator for the new Mat 1364659c6bb0SJed Brown . nr - number of nested row blocks 13650298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1366659c6bb0SJed Brown . nc - number of nested column blocks 13670298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 13680298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1369659c6bb0SJed Brown 1370659c6bb0SJed Brown Output Parameter: 1371659c6bb0SJed Brown . B - new matrix 1372659c6bb0SJed Brown 1373659c6bb0SJed Brown Level: advanced 1374659c6bb0SJed Brown 1375950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST 1376659c6bb0SJed Brown @*/ 13777087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1378d8588912SDave May { 1379d8588912SDave May Mat A; 1380d8588912SDave May PetscErrorCode ierr; 1381d8588912SDave May 1382d8588912SDave May PetscFunctionBegin; 1383c8883902SJed Brown *B = 0; 1384d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1385c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 13867ae8954aSSatish Balay ierr = MatSetUp(A);CHKERRQ(ierr); 1387c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1388d8588912SDave May *B = A; 1389d8588912SDave May PetscFunctionReturn(0); 1390d8588912SDave May } 1391659c6bb0SJed Brown 1392629c3df2SDmitry Karpeev #undef __FUNCT__ 1393629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ" 13948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1395629c3df2SDmitry Karpeev { 1396629c3df2SDmitry Karpeev PetscErrorCode ierr; 1397629c3df2SDmitry Karpeev Mat_Nest *nest = (Mat_Nest*)A->data; 1398629c3df2SDmitry Karpeev PetscInt m,n,M,N,i,j,k,*dnnz,*onnz; 1399629c3df2SDmitry Karpeev Mat C; 1400629c3df2SDmitry Karpeev 1401629c3df2SDmitry Karpeev PetscFunctionBegin; 1402629c3df2SDmitry Karpeev ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1403629c3df2SDmitry Karpeev ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1404629c3df2SDmitry Karpeev switch (reuse) { 1405629c3df2SDmitry Karpeev case MAT_INITIAL_MATRIX: 1406ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 1407629c3df2SDmitry Karpeev ierr = MatSetType(C,newtype);CHKERRQ(ierr); 1408629c3df2SDmitry Karpeev ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 1409629c3df2SDmitry Karpeev *newmat = C; 1410629c3df2SDmitry Karpeev break; 1411629c3df2SDmitry Karpeev case MAT_REUSE_MATRIX: 1412629c3df2SDmitry Karpeev C = *newmat; 1413629c3df2SDmitry Karpeev break; 1414ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse"); 1415629c3df2SDmitry Karpeev } 1416629c3df2SDmitry Karpeev 1417629c3df2SDmitry Karpeev /* Preallocation */ 1418785e854fSJed Brown ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr); 1419629c3df2SDmitry Karpeev onnz = dnnz + m; 1420629c3df2SDmitry Karpeev for (k=0; k<m; k++) { 1421629c3df2SDmitry Karpeev dnnz[k] = 0; 1422629c3df2SDmitry Karpeev onnz[k] = 0; 1423629c3df2SDmitry Karpeev } 1424629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1425629c3df2SDmitry Karpeev IS bNis; 1426629c3df2SDmitry Karpeev PetscInt bN; 1427629c3df2SDmitry Karpeev const PetscInt *bNindices; 1428629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1429629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1430629c3df2SDmitry Karpeev ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr); 1431629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1432629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1433629c3df2SDmitry Karpeev PetscSF bmsf; 1434629c3df2SDmitry Karpeev PetscSFNode *bmedges; 1435629c3df2SDmitry Karpeev Mat B; 143634986ce2SMatthew G Knepley PetscInt bm, *bmdnnz, br; 1437629c3df2SDmitry Karpeev const PetscInt *bmindices; 1438629c3df2SDmitry Karpeev B = nest->m[i][j]; 1439629c3df2SDmitry Karpeev if (!B) continue; 1440629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1441629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1442ce94432eSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr); 1443785e854fSJed Brown ierr = PetscMalloc1(2*bm,&bmedges);CHKERRQ(ierr); 1444785e854fSJed Brown ierr = PetscMalloc1(2*bm,&bmdnnz);CHKERRQ(ierr); 1445629c3df2SDmitry Karpeev for (k = 0; k < 2*bm; ++k) bmdnnz[k] = 0; 1446629c3df2SDmitry Karpeev /* 1447629c3df2SDmitry Karpeev Locate the owners for all of the locally-owned global row indices for this row block. 1448629c3df2SDmitry Karpeev These determine the roots of PetscSF used to communicate preallocation data to row owners. 1449629c3df2SDmitry Karpeev The roots correspond to the dnnz and onnz entries; thus, there are two roots per row. 1450629c3df2SDmitry Karpeev */ 1451629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1452a4b3d3acSMatthew G Knepley PetscInt row = bmindices[br], rowowner = 0, brncols, col, colowner = 0; 1453629c3df2SDmitry Karpeev const PetscInt *brcols; 1454a4b3d3acSMatthew G Knepley PetscInt rowrel = 0; /* row's relative index on its owner rank */ 1455629c3df2SDmitry Karpeev PetscInt rowownerm; /* local row size on row's owning rank. */ 145626fbe8dcSKarl Rupp 1457629c3df2SDmitry Karpeev ierr = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr); 1458629c3df2SDmitry Karpeev rowownerm = A->rmap->range[rowowner+1]-A->rmap->range[rowowner]; 145926fbe8dcSKarl Rupp 1460629c3df2SDmitry Karpeev bmedges[br].rank = rowowner; bmedges[br].index = rowrel; /* edge from bmdnnz to dnnz */ 1461629c3df2SDmitry Karpeev bmedges[br].rank = rowowner; bmedges[br].index = rowrel+rowownerm; /* edge from bmonnz to onnz */ 1462629c3df2SDmitry Karpeev /* Now actually compute the data -- bmdnnz and bmonnz by looking at the global columns in the br row of this block. */ 1463629c3df2SDmitry Karpeev /* Note that this is not a pessimistic bound only because we assume the index sets embedding the blocks do not overlap. */ 14640298fd71SBarry Smith ierr = MatGetRow(B,br,&brncols,&brcols,NULL);CHKERRQ(ierr); 1465629c3df2SDmitry Karpeev for (k=0; k<brncols; k++) { 1466629c3df2SDmitry Karpeev col = bNindices[brcols[k]]; 14670298fd71SBarry Smith ierr = PetscLayoutFindOwnerIndex(A->cmap,col,&colowner,NULL);CHKERRQ(ierr); 1468629c3df2SDmitry Karpeev if (colowner == rowowner) bmdnnz[br]++; 1469629c3df2SDmitry Karpeev else onnz[br]++; 1470629c3df2SDmitry Karpeev } 14710298fd71SBarry Smith ierr = MatRestoreRow(B,br,&brncols,&brcols,NULL);CHKERRQ(ierr); 1472629c3df2SDmitry Karpeev } 1473629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1474629c3df2SDmitry Karpeev /* bsf will have to take care of disposing of bedges. */ 14750298fd71SBarry Smith ierr = PetscSFSetGraph(bmsf,m,2*bm,NULL,PETSC_COPY_VALUES,bmedges,PETSC_OWN_POINTER);CHKERRQ(ierr); 1476629c3df2SDmitry Karpeev ierr = PetscSFReduceBegin(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr); 1477629c3df2SDmitry Karpeev ierr = PetscSFReduceEnd(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr); 1478629c3df2SDmitry Karpeev ierr = PetscFree(bmdnnz);CHKERRQ(ierr); 1479629c3df2SDmitry Karpeev ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr); 1480629c3df2SDmitry Karpeev } 148122d28d08SBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1482629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1483629c3df2SDmitry Karpeev } 1484629c3df2SDmitry Karpeev ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr); 1485629c3df2SDmitry Karpeev ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr); 1486629c3df2SDmitry Karpeev ierr = PetscFree(dnnz);CHKERRQ(ierr); 1487629c3df2SDmitry Karpeev 1488629c3df2SDmitry Karpeev /* Fill by row */ 1489629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1490629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1491629c3df2SDmitry Karpeev IS bNis; 1492629c3df2SDmitry Karpeev PetscInt bN; 1493629c3df2SDmitry Karpeev const PetscInt *bNindices; 1494629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1495629c3df2SDmitry Karpeev ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr); 1496629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1497629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1498629c3df2SDmitry Karpeev Mat B; 1499629c3df2SDmitry Karpeev PetscInt bm, br; 1500629c3df2SDmitry Karpeev const PetscInt *bmindices; 1501629c3df2SDmitry Karpeev B = nest->m[i][j]; 1502629c3df2SDmitry Karpeev if (!B) continue; 1503629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1504629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1505629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1506629c3df2SDmitry Karpeev PetscInt row = bmindices[br], brncols, *cols; 1507629c3df2SDmitry Karpeev const PetscInt *brcols; 1508629c3df2SDmitry Karpeev const PetscScalar *brcoldata; 1509629c3df2SDmitry Karpeev ierr = MatGetRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1510785e854fSJed Brown ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr); 151126fbe8dcSKarl Rupp for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]]; 1512629c3df2SDmitry Karpeev /* 1513629c3df2SDmitry Karpeev Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match. 1514629c3df2SDmitry Karpeev Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES. 1515629c3df2SDmitry Karpeev */ 1516a2ea699eSBarry Smith ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr); 1517629c3df2SDmitry Karpeev ierr = MatRestoreRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1518629c3df2SDmitry Karpeev ierr = PetscFree(cols);CHKERRQ(ierr); 1519629c3df2SDmitry Karpeev } 1520629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1521629c3df2SDmitry Karpeev } 1522a2ea699eSBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1523629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1524629c3df2SDmitry Karpeev } 1525629c3df2SDmitry Karpeev ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1526629c3df2SDmitry Karpeev ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1527629c3df2SDmitry Karpeev PetscFunctionReturn(0); 1528629c3df2SDmitry Karpeev } 1529629c3df2SDmitry Karpeev 1530659c6bb0SJed Brown /*MC 1531659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1532659c6bb0SJed Brown 1533659c6bb0SJed Brown Level: intermediate 1534659c6bb0SJed Brown 1535659c6bb0SJed Brown Notes: 1536659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1537659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1538950540a4SJed Brown It is usually used with DMComposite and DMCreateMatrix() 1539659c6bb0SJed Brown 1540659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1541659c6bb0SJed Brown M*/ 1542c8883902SJed Brown #undef __FUNCT__ 1543c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 15448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A) 1545c8883902SJed Brown { 1546c8883902SJed Brown Mat_Nest *s; 1547c8883902SJed Brown PetscErrorCode ierr; 1548c8883902SJed Brown 1549c8883902SJed Brown PetscFunctionBegin; 1550b00a9115SJed Brown ierr = PetscNewLog(A,&s);CHKERRQ(ierr); 1551c8883902SJed Brown A->data = (void*)s; 1552e7c19651SJed Brown 1553e7c19651SJed Brown s->nr = -1; 1554e7c19651SJed Brown s->nc = -1; 15550298fd71SBarry Smith s->m = NULL; 1556e7c19651SJed Brown s->splitassembly = PETSC_FALSE; 1557c8883902SJed Brown 1558c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 155926fbe8dcSKarl Rupp 1560c8883902SJed Brown A->ops->mult = MatMult_Nest; 15619194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 1562c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 15639194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1564c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1565c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1566c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1567*c222c20dSDavid Ham A->ops->copy = MatCopy_Nest; 1568c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1569c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1570c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1571c8883902SJed Brown A->ops->view = MatView_Nest; 1572c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1573c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1574c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1575429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 1576429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 1577a061e289SJed Brown A->ops->scale = MatScale_Nest; 1578a061e289SJed Brown A->ops->shift = MatShift_Nest; 1579c8883902SJed Brown 1580c8883902SJed Brown A->spptr = 0; 1581c8883902SJed Brown A->assembled = PETSC_FALSE; 1582c8883902SJed Brown 1583c8883902SJed Brown /* expose Nest api's */ 1584bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C", MatNestGetSubMat_Nest);CHKERRQ(ierr); 1585bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1586bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1587bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C", MatNestGetSize_Nest);CHKERRQ(ierr); 1588bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C", MatNestGetISs_Nest);CHKERRQ(ierr); 1589bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 1590bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C", MatNestSetVecType_Nest);CHKERRQ(ierr); 1591bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1592c8883902SJed Brown 1593c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1594c8883902SJed Brown PetscFunctionReturn(0); 1595c8883902SJed Brown } 1596