1d8588912SDave May 2ec9811eeSJed Brown #include "matnestimpl.h" /*I "petscmat.h" I*/ 3d8588912SDave May 4c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]); 57874fa86SDave May static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left); 6c8883902SJed Brown 7d8588912SDave May /* private functions */ 8d8588912SDave May #undef __FUNCT__ 98188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private" 108188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N) 11d8588912SDave May { 12d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 138188e55aSJed Brown PetscInt i,j; 14d8588912SDave May PetscErrorCode ierr; 15d8588912SDave May 16d8588912SDave May PetscFunctionBegin; 178188e55aSJed Brown *m = *n = *M = *N = 0; 188188e55aSJed Brown for (i=0; i<bA->nr; i++) { /* rows */ 198188e55aSJed Brown PetscInt sm,sM; 208188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr); 218188e55aSJed Brown ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr); 228188e55aSJed Brown *m += sm; 238188e55aSJed Brown *M += sM; 24d8588912SDave May } 258188e55aSJed Brown for (j=0; j<bA->nc; j++) { /* cols */ 268188e55aSJed Brown PetscInt sn,sN; 278188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr); 288188e55aSJed Brown ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr); 298188e55aSJed Brown *n += sn; 308188e55aSJed Brown *N += sN; 31d8588912SDave May } 32d8588912SDave May PetscFunctionReturn(0); 33d8588912SDave May } 34d8588912SDave May 35d8588912SDave May /* operations */ 36d8588912SDave May #undef __FUNCT__ 37d8588912SDave May #define __FUNCT__ "MatMult_Nest" 38207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 39d8588912SDave May { 40d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 41207556f9SJed Brown Vec *bx = bA->right,*by = bA->left; 42207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 43d8588912SDave May PetscErrorCode ierr; 44d8588912SDave May 45d8588912SDave May PetscFunctionBegin; 46207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 47207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 48207556f9SJed Brown for (i=0; i<nr; i++) { 49d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 50207556f9SJed Brown for (j=0; j<nc; j++) { 51207556f9SJed Brown if (!bA->m[i][j]) continue; 52d8588912SDave May /* y[i] <- y[i] + A[i][j] * x[j] */ 53d8588912SDave May ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 54d8588912SDave May } 55d8588912SDave May } 56207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 57207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 58d8588912SDave May PetscFunctionReturn(0); 59d8588912SDave May } 60d8588912SDave May 61d8588912SDave May #undef __FUNCT__ 629194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest" 639194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z) 649194d70fSJed Brown { 659194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 669194d70fSJed Brown Vec *bx = bA->right,*bz = bA->left; 679194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 689194d70fSJed Brown PetscErrorCode ierr; 699194d70fSJed Brown 709194d70fSJed Brown PetscFunctionBegin; 719194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 729194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 739194d70fSJed Brown for (i=0; i<nr; i++) { 749194d70fSJed Brown if (y != z) { 759194d70fSJed Brown Vec by; 769194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 779194d70fSJed Brown ierr = VecCopy(by,bz[i]);CHKERRQ(ierr); 789194d70fSJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 799194d70fSJed Brown } 809194d70fSJed Brown for (j=0; j<nc; j++) { 819194d70fSJed Brown if (!bA->m[i][j]) continue; 829194d70fSJed Brown /* y[i] <- y[i] + A[i][j] * x[j] */ 839194d70fSJed Brown ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr); 849194d70fSJed Brown } 859194d70fSJed Brown } 869194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 879194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 889194d70fSJed Brown PetscFunctionReturn(0); 899194d70fSJed Brown } 909194d70fSJed Brown 919194d70fSJed Brown #undef __FUNCT__ 92d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest" 93207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 94d8588912SDave May { 95d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 96207556f9SJed Brown Vec *bx = bA->left,*by = bA->right; 97207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 98d8588912SDave May PetscErrorCode ierr; 99d8588912SDave May 100d8588912SDave May PetscFunctionBegin; 101609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 102609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 103207556f9SJed Brown for (j=0; j<nc; j++) { 104609e31cbSJed Brown ierr = VecZeroEntries(by[j]);CHKERRQ(ierr); 105609e31cbSJed Brown for (i=0; i<nr; i++) { 106207556f9SJed Brown if (!bA->m[j][i]) continue; 107609e31cbSJed Brown /* y[j] <- y[j] + (A[i][j])^T * x[i] */ 108609e31cbSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr); 109d8588912SDave May } 110d8588912SDave May } 111609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 112609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 113d8588912SDave May PetscFunctionReturn(0); 114d8588912SDave May } 115d8588912SDave May 116d8588912SDave May #undef __FUNCT__ 1179194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest" 1189194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z) 1199194d70fSJed Brown { 1209194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 1219194d70fSJed Brown Vec *bx = bA->left,*bz = bA->right; 1229194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 1239194d70fSJed Brown PetscErrorCode ierr; 1249194d70fSJed Brown 1259194d70fSJed Brown PetscFunctionBegin; 1269194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 1279194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 1289194d70fSJed Brown for (j=0; j<nc; j++) { 1299194d70fSJed Brown if (y != z) { 1309194d70fSJed Brown Vec by; 1319194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 1329194d70fSJed Brown ierr = VecCopy(by,bz[j]);CHKERRQ(ierr); 1339194d70fSJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 1349194d70fSJed Brown } 1359194d70fSJed Brown for (i=0; i<nr; i++) { 1369194d70fSJed Brown if (!bA->m[j][i]) continue; 1379194d70fSJed Brown /* z[j] <- y[j] + (A[i][j])^T * x[i] */ 1389194d70fSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr); 1399194d70fSJed Brown } 1409194d70fSJed Brown } 1419194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 1429194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 1439194d70fSJed Brown PetscFunctionReturn(0); 1449194d70fSJed Brown } 1459194d70fSJed Brown 1469194d70fSJed Brown #undef __FUNCT__ 147e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList" 148e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 149e2d7f03fSJed Brown { 150e2d7f03fSJed Brown PetscErrorCode ierr; 151e2d7f03fSJed Brown IS *lst = *list; 152e2d7f03fSJed Brown PetscInt i; 153e2d7f03fSJed Brown 154e2d7f03fSJed Brown PetscFunctionBegin; 155e2d7f03fSJed Brown if (!lst) PetscFunctionReturn(0); 1566bf464f9SBarry Smith for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);} 157e2d7f03fSJed Brown ierr = PetscFree(lst);CHKERRQ(ierr); 158e2d7f03fSJed Brown *list = PETSC_NULL; 159e2d7f03fSJed Brown PetscFunctionReturn(0); 160e2d7f03fSJed Brown } 161e2d7f03fSJed Brown 162e2d7f03fSJed Brown #undef __FUNCT__ 163d8588912SDave May #define __FUNCT__ "MatDestroy_Nest" 164207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A) 165d8588912SDave May { 166d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 167d8588912SDave May PetscInt i,j; 168d8588912SDave May PetscErrorCode ierr; 169d8588912SDave May 170d8588912SDave May PetscFunctionBegin; 171d8588912SDave May /* release the matrices and the place holders */ 172e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 173e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 174e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 175e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 176d8588912SDave May 177d8588912SDave May ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 178d8588912SDave May ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 179d8588912SDave May 180207556f9SJed Brown ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 181207556f9SJed Brown 182d8588912SDave May /* release the matrices and the place holders */ 183d8588912SDave May if (vs->m) { 184d8588912SDave May for (i=0; i<vs->nr; i++) { 185d8588912SDave May for (j=0; j<vs->nc; j++) { 1866bf464f9SBarry Smith ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr); 187d8588912SDave May } 188d8588912SDave May ierr = PetscFree( vs->m[i] );CHKERRQ(ierr); 189d8588912SDave May } 190d8588912SDave May ierr = PetscFree(vs->m);CHKERRQ(ierr); 191d8588912SDave May } 192bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 193d8588912SDave May 194c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "",0);CHKERRQ(ierr); 1950782ca92SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "",0);CHKERRQ(ierr); 196c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "",0);CHKERRQ(ierr); 197c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "",0);CHKERRQ(ierr); 198c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "",0);CHKERRQ(ierr); 199c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "",0);CHKERRQ(ierr); 200d8588912SDave May PetscFunctionReturn(0); 201d8588912SDave May } 202d8588912SDave May 203d8588912SDave May #undef __FUNCT__ 204d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest" 205207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 206d8588912SDave May { 207d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 208d8588912SDave May PetscInt i,j; 209d8588912SDave May PetscErrorCode ierr; 210d8588912SDave May 211d8588912SDave May PetscFunctionBegin; 212d8588912SDave May for (i=0; i<vs->nr; i++) { 213d8588912SDave May for (j=0; j<vs->nc; j++) { 214d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); } 215d8588912SDave May } 216d8588912SDave May } 217d8588912SDave May PetscFunctionReturn(0); 218d8588912SDave May } 219d8588912SDave May 220d8588912SDave May #undef __FUNCT__ 221d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 222207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 223d8588912SDave May { 224d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 225d8588912SDave May PetscInt i,j; 226d8588912SDave May PetscErrorCode ierr; 227d8588912SDave May 228d8588912SDave May PetscFunctionBegin; 229d8588912SDave May for (i=0; i<vs->nr; i++) { 230d8588912SDave May for (j=0; j<vs->nc; j++) { 231d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); } 232d8588912SDave May } 233d8588912SDave May } 234d8588912SDave May PetscFunctionReturn(0); 235d8588912SDave May } 236d8588912SDave May 237d8588912SDave May #undef __FUNCT__ 238f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 239f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 240d8588912SDave May { 241207556f9SJed Brown PetscErrorCode ierr; 242f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 243f349c1fdSJed Brown PetscInt j; 244f349c1fdSJed Brown Mat sub; 245d8588912SDave May 246d8588912SDave May PetscFunctionBegin; 2478188e55aSJed Brown sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */ 248f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 2498188e55aSJed Brown if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 250f349c1fdSJed Brown *B = sub; 251f349c1fdSJed Brown PetscFunctionReturn(0); 252d8588912SDave May } 253d8588912SDave May 254f349c1fdSJed Brown #undef __FUNCT__ 255f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 256f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 257f349c1fdSJed Brown { 258207556f9SJed Brown PetscErrorCode ierr; 259f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 260f349c1fdSJed Brown PetscInt i; 261f349c1fdSJed Brown Mat sub; 262f349c1fdSJed Brown 263f349c1fdSJed Brown PetscFunctionBegin; 2648188e55aSJed Brown sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */ 265f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 2668188e55aSJed Brown if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 267f349c1fdSJed Brown *B = sub; 268f349c1fdSJed Brown PetscFunctionReturn(0); 269d8588912SDave May } 270d8588912SDave May 271f349c1fdSJed Brown #undef __FUNCT__ 272f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 273f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 274f349c1fdSJed Brown { 275f349c1fdSJed Brown PetscErrorCode ierr; 276f349c1fdSJed Brown PetscInt i; 277f349c1fdSJed Brown PetscBool flg; 278f349c1fdSJed Brown 279f349c1fdSJed Brown PetscFunctionBegin; 280f349c1fdSJed Brown PetscValidPointer(list,3); 281f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 282f349c1fdSJed Brown PetscValidIntPointer(found,5); 283f349c1fdSJed Brown *found = -1; 284f349c1fdSJed Brown for (i=0; i<n; i++) { 285207556f9SJed Brown if (!list[i]) continue; 286f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 287f349c1fdSJed Brown if (flg) { 288f349c1fdSJed Brown *found = i; 289f349c1fdSJed Brown PetscFunctionReturn(0); 290f349c1fdSJed Brown } 291f349c1fdSJed Brown } 292f349c1fdSJed Brown SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set"); 293f349c1fdSJed Brown PetscFunctionReturn(0); 294f349c1fdSJed Brown } 295f349c1fdSJed Brown 296f349c1fdSJed Brown #undef __FUNCT__ 2978188e55aSJed Brown #define __FUNCT__ "MatNestGetRow" 2988188e55aSJed Brown /* Get a block row as a new MatNest */ 2998188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 3008188e55aSJed Brown { 3018188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3028188e55aSJed Brown char keyname[256]; 3038188e55aSJed Brown PetscErrorCode ierr; 3048188e55aSJed Brown 3058188e55aSJed Brown PetscFunctionBegin; 3068188e55aSJed Brown *B = PETSC_NULL; 3078188e55aSJed Brown ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr); 3088188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 3098188e55aSJed Brown if (*B) PetscFunctionReturn(0); 3108188e55aSJed Brown 3118188e55aSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 3128188e55aSJed Brown (*B)->assembled = A->assembled; 3138188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 3148188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 3158188e55aSJed Brown PetscFunctionReturn(0); 3168188e55aSJed Brown } 3178188e55aSJed Brown 3188188e55aSJed Brown #undef __FUNCT__ 319f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 320f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 321f349c1fdSJed Brown { 322f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3238188e55aSJed Brown PetscErrorCode ierr; 3246b3a5b13SJed Brown PetscInt row,col; 3258188e55aSJed Brown PetscBool same,isFullCol; 326f349c1fdSJed Brown 327f349c1fdSJed Brown PetscFunctionBegin; 3288188e55aSJed Brown /* Check if full column space. This is a hack */ 3298188e55aSJed Brown isFullCol = PETSC_FALSE; 3308188e55aSJed Brown ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 3318188e55aSJed Brown if (same) { 3328188e55aSJed Brown PetscInt n,first,step; 3338188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 3348188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 3358188e55aSJed Brown if (A->cmap->n == n && A->cmap->rstart == first && step == 1) isFullCol = PETSC_TRUE; 3368188e55aSJed Brown } 3378188e55aSJed Brown 3388188e55aSJed Brown if (isFullCol) { 3398188e55aSJed Brown PetscInt row; 3408188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 3418188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 3428188e55aSJed Brown } else { 343f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 344f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 345f349c1fdSJed Brown *B = vs->m[row][col]; 3468188e55aSJed Brown } 347f349c1fdSJed Brown PetscFunctionReturn(0); 348f349c1fdSJed Brown } 349f349c1fdSJed Brown 350f349c1fdSJed Brown #undef __FUNCT__ 351f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 352f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 353f349c1fdSJed Brown { 354f349c1fdSJed Brown PetscErrorCode ierr; 355f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 356f349c1fdSJed Brown Mat sub; 357f349c1fdSJed Brown 358f349c1fdSJed Brown PetscFunctionBegin; 359f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 360f349c1fdSJed Brown switch (reuse) { 361f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 3627874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 363f349c1fdSJed Brown *B = sub; 364f349c1fdSJed Brown break; 365f349c1fdSJed Brown case MAT_REUSE_MATRIX: 366f349c1fdSJed Brown if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 367f349c1fdSJed Brown break; 368f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 369f349c1fdSJed Brown break; 370f349c1fdSJed Brown } 371f349c1fdSJed Brown PetscFunctionReturn(0); 372f349c1fdSJed Brown } 373f349c1fdSJed Brown 374f349c1fdSJed Brown #undef __FUNCT__ 375f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 376f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 377f349c1fdSJed Brown { 378f349c1fdSJed Brown PetscErrorCode ierr; 379f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 380f349c1fdSJed Brown Mat sub; 381f349c1fdSJed Brown 382f349c1fdSJed Brown PetscFunctionBegin; 383f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 384f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 385f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 386f349c1fdSJed Brown *B = sub; 387d8588912SDave May PetscFunctionReturn(0); 388d8588912SDave May } 389d8588912SDave May 390d8588912SDave May #undef __FUNCT__ 391d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 392207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 393d8588912SDave May { 394d8588912SDave May PetscErrorCode ierr; 395f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 396f349c1fdSJed Brown Mat sub; 397d8588912SDave May 398d8588912SDave May PetscFunctionBegin; 399f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 400f349c1fdSJed Brown if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 401f349c1fdSJed Brown if (sub) { 402f349c1fdSJed Brown if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4036bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 404d8588912SDave May } 405d8588912SDave May PetscFunctionReturn(0); 406d8588912SDave May } 407d8588912SDave May 408d8588912SDave May #undef __FUNCT__ 4097874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 4107874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 4117874fa86SDave May { 4127874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4137874fa86SDave May PetscInt i; 4147874fa86SDave May PetscErrorCode ierr; 4157874fa86SDave May 4167874fa86SDave May PetscFunctionBegin; 4177874fa86SDave May for (i=0; i<bA->nr; i++) { 418429bac76SJed Brown Vec bv; 419429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4207874fa86SDave May if (bA->m[i][i]) { 421429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 4227874fa86SDave May } else { 423429bac76SJed Brown ierr = VecSet(bv,1.0);CHKERRQ(ierr); 4247874fa86SDave May } 425429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4267874fa86SDave May } 4277874fa86SDave May PetscFunctionReturn(0); 4287874fa86SDave May } 4297874fa86SDave May 4307874fa86SDave May #undef __FUNCT__ 4317874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 4327874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 4337874fa86SDave May { 4347874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 435429bac76SJed Brown Vec bl,*br; 4367874fa86SDave May PetscInt i,j; 4377874fa86SDave May PetscErrorCode ierr; 4387874fa86SDave May 4397874fa86SDave May PetscFunctionBegin; 440429bac76SJed Brown ierr = PetscMalloc(bA->nc*sizeof(Vec),&br);CHKERRQ(ierr); 441429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 4427874fa86SDave May for (i=0; i<bA->nr; i++) { 443429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4447874fa86SDave May for (j=0; j<bA->nc; j++) { 4457874fa86SDave May if (bA->m[i][j]) { 446429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 4477874fa86SDave May } 4487874fa86SDave May } 449*a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4507874fa86SDave May } 451429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 452429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 4537874fa86SDave May PetscFunctionReturn(0); 4547874fa86SDave May } 4557874fa86SDave May 4567874fa86SDave May #undef __FUNCT__ 457*a061e289SJed Brown #define __FUNCT__ "MatScale_Nest" 458*a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 459*a061e289SJed Brown { 460*a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 461*a061e289SJed Brown PetscInt i,j; 462*a061e289SJed Brown PetscErrorCode ierr; 463*a061e289SJed Brown 464*a061e289SJed Brown PetscFunctionBegin; 465*a061e289SJed Brown for (i=0; i<bA->nr; i++) { 466*a061e289SJed Brown for (j=0; j<bA->nc; j++) { 467*a061e289SJed Brown if (bA->m[i][j]) { 468*a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 469*a061e289SJed Brown } 470*a061e289SJed Brown } 471*a061e289SJed Brown } 472*a061e289SJed Brown PetscFunctionReturn(0); 473*a061e289SJed Brown } 474*a061e289SJed Brown 475*a061e289SJed Brown #undef __FUNCT__ 476*a061e289SJed Brown #define __FUNCT__ "MatShift_Nest" 477*a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 478*a061e289SJed Brown { 479*a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 480*a061e289SJed Brown PetscInt i; 481*a061e289SJed Brown PetscErrorCode ierr; 482*a061e289SJed Brown 483*a061e289SJed Brown PetscFunctionBegin; 484*a061e289SJed Brown for (i=0; i<bA->nr; i++) { 485*a061e289SJed Brown if (!bA->m[i][i]) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i); 486*a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 487*a061e289SJed Brown } 488*a061e289SJed Brown PetscFunctionReturn(0); 489*a061e289SJed Brown } 490*a061e289SJed Brown 491*a061e289SJed Brown #undef __FUNCT__ 492d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest" 493207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 494d8588912SDave May { 495d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 496d8588912SDave May Vec *L,*R; 497d8588912SDave May MPI_Comm comm; 498d8588912SDave May PetscInt i,j; 499d8588912SDave May PetscErrorCode ierr; 500d8588912SDave May 501d8588912SDave May PetscFunctionBegin; 502d8588912SDave May comm = ((PetscObject)A)->comm; 503d8588912SDave May if (right) { 504d8588912SDave May /* allocate R */ 505d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr); 506d8588912SDave May /* Create the right vectors */ 507d8588912SDave May for (j=0; j<bA->nc; j++) { 508d8588912SDave May for (i=0; i<bA->nr; i++) { 509d8588912SDave May if (bA->m[i][j]) { 510d8588912SDave May ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 511d8588912SDave May break; 512d8588912SDave May } 513d8588912SDave May } 514d8588912SDave May if (i==bA->nr) { 515d8588912SDave May /* have an empty column */ 516d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 517d8588912SDave May } 518d8588912SDave May } 519f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 520d8588912SDave May /* hand back control to the nest vector */ 521d8588912SDave May for (j=0; j<bA->nc; j++) { 5226bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 523d8588912SDave May } 524d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 525d8588912SDave May } 526d8588912SDave May 527d8588912SDave May if (left) { 528d8588912SDave May /* allocate L */ 529d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr); 530d8588912SDave May /* Create the left vectors */ 531d8588912SDave May for (i=0; i<bA->nr; i++) { 532d8588912SDave May for (j=0; j<bA->nc; j++) { 533d8588912SDave May if (bA->m[i][j]) { 534d8588912SDave May ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 535d8588912SDave May break; 536d8588912SDave May } 537d8588912SDave May } 538d8588912SDave May if (j==bA->nc) { 539d8588912SDave May /* have an empty row */ 540d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 541d8588912SDave May } 542d8588912SDave May } 543d8588912SDave May 544f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 545d8588912SDave May for (i=0; i<bA->nr; i++) { 5466bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 547d8588912SDave May } 548d8588912SDave May 549d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 550d8588912SDave May } 551d8588912SDave May PetscFunctionReturn(0); 552d8588912SDave May } 553d8588912SDave May 554d8588912SDave May #undef __FUNCT__ 555d8588912SDave May #define __FUNCT__ "MatView_Nest" 556207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 557d8588912SDave May { 558d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 559d8588912SDave May PetscBool isascii; 560d8588912SDave May PetscInt i,j; 561d8588912SDave May PetscErrorCode ierr; 562d8588912SDave May 563d8588912SDave May PetscFunctionBegin; 564d8588912SDave May ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 565d8588912SDave May if (isascii) { 566d8588912SDave May 567d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n" ); 568d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push0 */ 569d8588912SDave May PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 570d8588912SDave May 571d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" ); 572d8588912SDave May for (i=0; i<bA->nr; i++) { 573d8588912SDave May for (j=0; j<bA->nc; j++) { 574d8588912SDave May const MatType type; 575d8588912SDave May const char *name; 576d8588912SDave May PetscInt NR,NC; 577d8588912SDave May PetscBool isNest = PETSC_FALSE; 578d8588912SDave May 579d8588912SDave May if (!bA->m[i][j]) { 580d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j); 581d8588912SDave May continue; 582d8588912SDave May } 583d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 584d8588912SDave May ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr); 585d8588912SDave May name = ((PetscObject)bA->m[i][j])->prefix; 586d8588912SDave May ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 587d8588912SDave May 588d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC); 589d8588912SDave May 590d8588912SDave May if (isNest) { 591d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push1 */ 592d8588912SDave May ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr); 593d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop1 */ 594d8588912SDave May } 595d8588912SDave May } 596d8588912SDave May } 597d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 598d8588912SDave May } 599d8588912SDave May PetscFunctionReturn(0); 600d8588912SDave May } 601d8588912SDave May 602d8588912SDave May #undef __FUNCT__ 603d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 604207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 605d8588912SDave May { 606d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 607d8588912SDave May PetscInt i,j; 608d8588912SDave May PetscErrorCode ierr; 609d8588912SDave May 610d8588912SDave May PetscFunctionBegin; 611d8588912SDave May for (i=0; i<bA->nr; i++) { 612d8588912SDave May for (j=0; j<bA->nc; j++) { 613d8588912SDave May if (!bA->m[i][j]) continue; 614d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 615d8588912SDave May } 616d8588912SDave May } 617d8588912SDave May PetscFunctionReturn(0); 618d8588912SDave May } 619d8588912SDave May 620d8588912SDave May #undef __FUNCT__ 621d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 622207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 623d8588912SDave May { 624d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 625841e96a3SJed Brown Mat *b; 626841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 627d8588912SDave May PetscErrorCode ierr; 628d8588912SDave May 629d8588912SDave May PetscFunctionBegin; 630841e96a3SJed Brown ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 631841e96a3SJed Brown for (i=0; i<nr; i++) { 632841e96a3SJed Brown for (j=0; j<nc; j++) { 633841e96a3SJed Brown if (bA->m[i][j]) { 634841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 635841e96a3SJed Brown } else { 636841e96a3SJed Brown b[i*nc+j] = PETSC_NULL; 637d8588912SDave May } 638d8588912SDave May } 639d8588912SDave May } 640f349c1fdSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 641841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 642841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 6436bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 644d8588912SDave May } 645d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 646d8588912SDave May 647841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 648841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 649d8588912SDave May PetscFunctionReturn(0); 650d8588912SDave May } 651d8588912SDave May 652d8588912SDave May /* nest api */ 653d8588912SDave May EXTERN_C_BEGIN 654d8588912SDave May #undef __FUNCT__ 655d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 656d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 657d8588912SDave May { 658d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 659d8588912SDave May PetscFunctionBegin; 660d8588912SDave May if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 661d8588912SDave May if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 662d8588912SDave May *mat = bA->m[idxm][jdxm]; 663d8588912SDave May PetscFunctionReturn(0); 664d8588912SDave May } 665d8588912SDave May EXTERN_C_END 666d8588912SDave May 667d8588912SDave May #undef __FUNCT__ 668d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 669d8588912SDave May /*@C 670d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 671d8588912SDave May 672d8588912SDave May Not collective 673d8588912SDave May 674d8588912SDave May Input Parameters: 675629881c0SJed Brown + A - nest matrix 676d8588912SDave May . idxm - index of the matrix within the nest matrix 677629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 678d8588912SDave May 679d8588912SDave May Output Parameter: 680d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 681d8588912SDave May 682d8588912SDave May Level: developer 683d8588912SDave May 684d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 685d8588912SDave May @*/ 6867087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 687d8588912SDave May { 688699a902aSJed Brown PetscErrorCode ierr; 689d8588912SDave May 690d8588912SDave May PetscFunctionBegin; 691699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 692d8588912SDave May PetscFunctionReturn(0); 693d8588912SDave May } 694d8588912SDave May 695d8588912SDave May EXTERN_C_BEGIN 696d8588912SDave May #undef __FUNCT__ 6970782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 6980782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 6990782ca92SJed Brown { 7000782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 7010782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 7020782ca92SJed Brown PetscErrorCode ierr; 7030782ca92SJed Brown 7040782ca92SJed Brown PetscFunctionBegin; 7050782ca92SJed Brown if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 7060782ca92SJed Brown if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 7070782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 7080782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 7090782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 7100782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 7110782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 7120782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 7130782ca92SJed Brown if (M != Mi || N != Ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni); 7140782ca92SJed Brown if (m != mi || n != ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni); 7150782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 7160782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 7170782ca92SJed Brown bA->m[idxm][jdxm] = mat; 7180782ca92SJed Brown PetscFunctionReturn(0); 7190782ca92SJed Brown } 7200782ca92SJed Brown EXTERN_C_END 7210782ca92SJed Brown 7220782ca92SJed Brown #undef __FUNCT__ 7230782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 7240782ca92SJed Brown /*@C 7250782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 7260782ca92SJed Brown 7270782ca92SJed Brown Logically collective on the submatrix communicator 7280782ca92SJed Brown 7290782ca92SJed Brown Input Parameters: 7300782ca92SJed Brown + A - nest matrix 7310782ca92SJed Brown . idxm - index of the matrix within the nest matrix 7320782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 7330782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 7340782ca92SJed Brown 7350782ca92SJed Brown Notes: 7360782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 7370782ca92SJed Brown 7380782ca92SJed Brown This increments the reference count of the submatrix. 7390782ca92SJed Brown 7400782ca92SJed Brown Level: developer 7410782ca92SJed Brown 7420782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 7430782ca92SJed Brown @*/ 7440782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 7450782ca92SJed Brown { 7460782ca92SJed Brown PetscErrorCode ierr; 7470782ca92SJed Brown 7480782ca92SJed Brown PetscFunctionBegin; 7490782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 7500782ca92SJed Brown PetscFunctionReturn(0); 7510782ca92SJed Brown } 7520782ca92SJed Brown 7530782ca92SJed Brown EXTERN_C_BEGIN 7540782ca92SJed Brown #undef __FUNCT__ 755d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 756d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 757d8588912SDave May { 758d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 759d8588912SDave May PetscFunctionBegin; 760d8588912SDave May if (M) { *M = bA->nr; } 761d8588912SDave May if (N) { *N = bA->nc; } 762d8588912SDave May if (mat) { *mat = bA->m; } 763d8588912SDave May PetscFunctionReturn(0); 764d8588912SDave May } 765d8588912SDave May EXTERN_C_END 766d8588912SDave May 767d8588912SDave May #undef __FUNCT__ 768d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 769d8588912SDave May /*@C 770d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 771d8588912SDave May 772d8588912SDave May Not collective 773d8588912SDave May 774d8588912SDave May Input Parameters: 775629881c0SJed Brown . A - nest matrix 776d8588912SDave May 777d8588912SDave May Output Parameter: 778629881c0SJed Brown + M - number of rows in the nest matrix 779d8588912SDave May . N - number of cols in the nest matrix 780629881c0SJed Brown - mat - 2d array of matrices 781d8588912SDave May 782d8588912SDave May Notes: 783d8588912SDave May 784d8588912SDave May The user should not free the array mat. 785d8588912SDave May 786d8588912SDave May Level: developer 787d8588912SDave May 788d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 789d8588912SDave May @*/ 7907087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 791d8588912SDave May { 792699a902aSJed Brown PetscErrorCode ierr; 793d8588912SDave May 794d8588912SDave May PetscFunctionBegin; 795699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 796d8588912SDave May PetscFunctionReturn(0); 797d8588912SDave May } 798d8588912SDave May 799d8588912SDave May EXTERN_C_BEGIN 800d8588912SDave May #undef __FUNCT__ 801d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 8027087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 803d8588912SDave May { 804d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 805d8588912SDave May 806d8588912SDave May PetscFunctionBegin; 807d8588912SDave May if (M) { *M = bA->nr; } 808d8588912SDave May if (N) { *N = bA->nc; } 809d8588912SDave May PetscFunctionReturn(0); 810d8588912SDave May } 811d8588912SDave May EXTERN_C_END 812d8588912SDave May 813d8588912SDave May #undef __FUNCT__ 814d8588912SDave May #define __FUNCT__ "MatNestGetSize" 815d8588912SDave May /*@C 816d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 817d8588912SDave May 818d8588912SDave May Not collective 819d8588912SDave May 820d8588912SDave May Input Parameters: 821d8588912SDave May . A - nest matrix 822d8588912SDave May 823d8588912SDave May Output Parameter: 824629881c0SJed Brown + M - number of rows in the nested mat 825629881c0SJed Brown - N - number of cols in the nested mat 826d8588912SDave May 827d8588912SDave May Notes: 828d8588912SDave May 829d8588912SDave May Level: developer 830d8588912SDave May 831d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 832d8588912SDave May @*/ 8337087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 834d8588912SDave May { 835699a902aSJed Brown PetscErrorCode ierr; 836d8588912SDave May 837d8588912SDave May PetscFunctionBegin; 838699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 839d8588912SDave May PetscFunctionReturn(0); 840d8588912SDave May } 841d8588912SDave May 842207556f9SJed Brown EXTERN_C_BEGIN 843207556f9SJed Brown #undef __FUNCT__ 844207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 8457087cfbeSBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,const VecType vtype) 846207556f9SJed Brown { 847207556f9SJed Brown PetscErrorCode ierr; 848207556f9SJed Brown PetscBool flg; 849207556f9SJed Brown 850207556f9SJed Brown PetscFunctionBegin; 851207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 852207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 853207556f9SJed Brown A->ops->getvecs = flg ? MatGetVecs_Nest : 0; 854207556f9SJed Brown PetscFunctionReturn(0); 855207556f9SJed Brown } 856207556f9SJed Brown EXTERN_C_END 857207556f9SJed Brown 858207556f9SJed Brown #undef __FUNCT__ 859207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 860207556f9SJed Brown /*@C 861207556f9SJed Brown MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 862207556f9SJed Brown 863207556f9SJed Brown Not collective 864207556f9SJed Brown 865207556f9SJed Brown Input Parameters: 866207556f9SJed Brown + A - nest matrix 867207556f9SJed Brown - vtype - type to use for creating vectors 868207556f9SJed Brown 869207556f9SJed Brown Notes: 870207556f9SJed Brown 871207556f9SJed Brown Level: developer 872207556f9SJed Brown 873207556f9SJed Brown .seealso: MatGetVecs() 874207556f9SJed Brown @*/ 8757087cfbeSBarry Smith PetscErrorCode MatNestSetVecType(Mat A,const VecType vtype) 876207556f9SJed Brown { 877207556f9SJed Brown PetscErrorCode ierr; 878207556f9SJed Brown 879207556f9SJed Brown PetscFunctionBegin; 880207556f9SJed Brown ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr); 881207556f9SJed Brown PetscFunctionReturn(0); 882207556f9SJed Brown } 883207556f9SJed Brown 884c8883902SJed Brown EXTERN_C_BEGIN 885d8588912SDave May #undef __FUNCT__ 886c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 887c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 888d8588912SDave May { 889c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 890c8883902SJed Brown PetscInt i,j,m,n,M,N; 891d8588912SDave May PetscErrorCode ierr; 892d8588912SDave May 893d8588912SDave May PetscFunctionBegin; 894c8883902SJed Brown s->nr = nr; 895c8883902SJed Brown s->nc = nc; 896d8588912SDave May 897c8883902SJed Brown /* Create space for submatrices */ 898c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 899c8883902SJed Brown for (i=0; i<nr; i++) { 900c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 901d8588912SDave May } 902c8883902SJed Brown for (i=0; i<nr; i++) { 903c8883902SJed Brown for (j=0; j<nc; j++) { 904c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 905c8883902SJed Brown if (a[i*nc+j]) { 906c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 907d8588912SDave May } 908d8588912SDave May } 909d8588912SDave May } 910d8588912SDave May 9118188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 912d8588912SDave May 913c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 914c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 915c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 916c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 917d8588912SDave May 9188188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 919d8588912SDave May 920c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 921c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 922c8883902SJed Brown ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 923c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 924c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 925c8883902SJed Brown ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 926c8883902SJed Brown 927c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 928c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 929c8883902SJed Brown 930c8883902SJed Brown ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 931c8883902SJed Brown ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr); 932c8883902SJed Brown ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr); 933d8588912SDave May PetscFunctionReturn(0); 934d8588912SDave May } 935c8883902SJed Brown EXTERN_C_END 936d8588912SDave May 937c8883902SJed Brown #undef __FUNCT__ 938c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 939c8883902SJed Brown /*@ 940c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 941c8883902SJed Brown 942c8883902SJed Brown Collective on Mat 943c8883902SJed Brown 944c8883902SJed Brown Input Parameter: 945c8883902SJed Brown + N - nested matrix 946c8883902SJed Brown . nr - number of nested row blocks 947c8883902SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 948c8883902SJed Brown . nc - number of nested column blocks 949c8883902SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 950c8883902SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 951c8883902SJed Brown 952c8883902SJed Brown Level: advanced 953c8883902SJed Brown 954c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 955c8883902SJed Brown @*/ 956c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 957c8883902SJed Brown { 958c8883902SJed Brown PetscErrorCode ierr; 959c8883902SJed Brown PetscInt i; 960c8883902SJed Brown 961c8883902SJed Brown PetscFunctionBegin; 962c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 963c8883902SJed Brown if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 964c8883902SJed Brown if (nr && is_row) { 965c8883902SJed Brown PetscValidPointer(is_row,3); 966c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 967c8883902SJed Brown } 968c8883902SJed Brown if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 9691664e352SJed Brown if (nc && is_col) { 970c8883902SJed Brown PetscValidPointer(is_col,5); 971c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 972c8883902SJed Brown } 973c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 974c8883902SJed 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); 975c8883902SJed Brown PetscFunctionReturn(0); 976c8883902SJed Brown } 977d8588912SDave May 978d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 979d8588912SDave May /* 980d8588912SDave May nprocessors = NP 981d8588912SDave May Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) ) 982d8588912SDave May proc 0: => (g_0,h_0,) 983d8588912SDave May proc 1: => (g_1,h_1,) 984d8588912SDave May ... 985d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 986d8588912SDave May 987d8588912SDave May proc 0: proc 1: proc nprocs-1: 988d8588912SDave May is[0] = ( 0,1,2,...,nlocal(g_0)-1 ) ( 0,1,...,nlocal(g_1)-1 ) ( 0,1,...,nlocal(g_NP-1) ) 989d8588912SDave May 990d8588912SDave May proc 0: 991d8588912SDave May is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 ) 992d8588912SDave May proc 1: 993d8588912SDave May is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 ) 994d8588912SDave May 995d8588912SDave May proc NP-1: 996d8588912SDave May is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 ) 997d8588912SDave May */ 998d8588912SDave May #undef __FUNCT__ 999d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1000841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1001d8588912SDave May { 1002e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 10038188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1004d8588912SDave May PetscErrorCode ierr; 10052ae74bdbSJed Brown Mat sub; 1006d8588912SDave May 1007d8588912SDave May PetscFunctionBegin; 10088188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 10098188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1010d8588912SDave May if (is_row) { /* valid IS is passed in */ 1011d8588912SDave May /* refs on is[] are incremeneted */ 1012e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1013d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1014e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1015d8588912SDave May } 10162ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 10178188e55aSJed Brown nsum = 0; 10188188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 10198188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 10208188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 10218188e55aSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 10228188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 10238188e55aSJed Brown nsum += n; 10248188e55aSJed Brown } 10258188e55aSJed Brown offset = 0; 10268188e55aSJed Brown ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1027e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1028f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 10292ae74bdbSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 10302ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1031e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1032e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 10332ae74bdbSJed Brown offset += n; 1034d8588912SDave May } 1035d8588912SDave May } 1036d8588912SDave May 1037d8588912SDave May if (is_col) { /* valid IS is passed in */ 1038d8588912SDave May /* refs on is[] are incremeneted */ 1039e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1040d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1041e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1042d8588912SDave May } 10432ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 10442ae74bdbSJed Brown offset = A->cmap->rstart; 10458188e55aSJed Brown nsum = 0; 10468188e55aSJed Brown for (j=0; j<vs->nc; j++) { 10478188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 10488188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 10498188e55aSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 10508188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 10518188e55aSJed Brown nsum += n; 10528188e55aSJed Brown } 10538188e55aSJed Brown offset = 0; 10548188e55aSJed Brown ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1055e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1056f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 10572ae74bdbSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 10582ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1059e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1060e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 10612ae74bdbSJed Brown offset += n; 1062d8588912SDave May } 1063d8588912SDave May } 1064e2d7f03fSJed Brown 1065e2d7f03fSJed Brown /* Set up the local ISs */ 1066e2d7f03fSJed Brown ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1067e2d7f03fSJed Brown ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1068e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1069e2d7f03fSJed Brown IS isloc; 10708188e55aSJed Brown ISLocalToGlobalMapping rmap = PETSC_NULL; 1071e2d7f03fSJed Brown PetscInt nlocal,bs; 1072e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 10738188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);} 1074207556f9SJed Brown if (rmap) { 1075e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1076e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1077e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1078e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1079207556f9SJed Brown } else { 1080207556f9SJed Brown nlocal = 0; 1081207556f9SJed Brown isloc = PETSC_NULL; 1082207556f9SJed Brown } 1083e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1084e2d7f03fSJed Brown offset += nlocal; 1085e2d7f03fSJed Brown } 10868188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1087e2d7f03fSJed Brown IS isloc; 10888188e55aSJed Brown ISLocalToGlobalMapping cmap = PETSC_NULL; 1089e2d7f03fSJed Brown PetscInt nlocal,bs; 1090e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 10918188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);} 1092207556f9SJed Brown if (cmap) { 1093e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1094e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1095e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1096e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1097207556f9SJed Brown } else { 1098207556f9SJed Brown nlocal = 0; 1099207556f9SJed Brown isloc = PETSC_NULL; 1100207556f9SJed Brown } 1101e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1102e2d7f03fSJed Brown offset += nlocal; 1103e2d7f03fSJed Brown } 11040189643fSJed Brown 11050189643fSJed Brown #if defined(PETSC_USE_DEBUG) 11060189643fSJed Brown for (i=0; i<vs->nr; i++) { 11070189643fSJed Brown for (j=0; j<vs->nc; j++) { 11080189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 11090189643fSJed Brown Mat B = vs->m[i][j]; 11100189643fSJed Brown if (!B) continue; 11110189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 11120189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 11130189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 11140189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 11150189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 11160189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 11170189643fSJed Brown if (M != Mi || N != Ni) SETERRQ6(((PetscObject)sub)->comm,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); 11180189643fSJed Brown if (m != mi || n != ni) SETERRQ6(((PetscObject)sub)->comm,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); 11190189643fSJed Brown } 11200189643fSJed Brown } 11210189643fSJed Brown #endif 1122*a061e289SJed Brown 1123*a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1124*a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1125*a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1126*a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1127*a061e289SJed Brown } 1128*a061e289SJed Brown } 1129*a061e289SJed Brown A->assembled = PETSC_TRUE; 1130d8588912SDave May PetscFunctionReturn(0); 1131d8588912SDave May } 1132d8588912SDave May 1133d8588912SDave May #undef __FUNCT__ 1134d8588912SDave May #define __FUNCT__ "MatCreateNest" 1135659c6bb0SJed Brown /*@ 1136659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1137659c6bb0SJed Brown 1138659c6bb0SJed Brown Collective on Mat 1139659c6bb0SJed Brown 1140659c6bb0SJed Brown Input Parameter: 1141659c6bb0SJed Brown + comm - Communicator for the new Mat 1142659c6bb0SJed Brown . nr - number of nested row blocks 1143659c6bb0SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1144659c6bb0SJed Brown . nc - number of nested column blocks 1145659c6bb0SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1146659c6bb0SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1147659c6bb0SJed Brown 1148659c6bb0SJed Brown Output Parameter: 1149659c6bb0SJed Brown . B - new matrix 1150659c6bb0SJed Brown 1151659c6bb0SJed Brown Level: advanced 1152659c6bb0SJed Brown 1153659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST 1154659c6bb0SJed Brown @*/ 11557087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1156d8588912SDave May { 1157d8588912SDave May Mat A; 1158d8588912SDave May PetscErrorCode ierr; 1159d8588912SDave May 1160d8588912SDave May PetscFunctionBegin; 1161c8883902SJed Brown *B = 0; 1162d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1163c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 1164c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1165d8588912SDave May *B = A; 1166d8588912SDave May PetscFunctionReturn(0); 1167d8588912SDave May } 1168659c6bb0SJed Brown 1169659c6bb0SJed Brown /*MC 1170659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1171659c6bb0SJed Brown 1172659c6bb0SJed Brown Level: intermediate 1173659c6bb0SJed Brown 1174659c6bb0SJed Brown Notes: 1175659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1176659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1177659c6bb0SJed Brown It is usually used with DMComposite and DMGetMatrix() 1178659c6bb0SJed Brown 1179659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1180659c6bb0SJed Brown M*/ 1181c8883902SJed Brown EXTERN_C_BEGIN 1182c8883902SJed Brown #undef __FUNCT__ 1183c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 1184c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A) 1185c8883902SJed Brown { 1186c8883902SJed Brown Mat_Nest *s; 1187c8883902SJed Brown PetscErrorCode ierr; 1188c8883902SJed Brown 1189c8883902SJed Brown PetscFunctionBegin; 1190c8883902SJed Brown ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr); 1191c8883902SJed Brown A->data = (void*)s; 1192c8883902SJed Brown s->nr = s->nc = -1; 1193c8883902SJed Brown s->m = PETSC_NULL; 1194c8883902SJed Brown 1195c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 1196c8883902SJed Brown A->ops->mult = MatMult_Nest; 11979194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 1198c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 11999194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1200c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1201c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1202c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1203c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1204c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1205c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1206c8883902SJed Brown A->ops->view = MatView_Nest; 1207c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1208c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1209c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1210429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 1211429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 1212*a061e289SJed Brown A->ops->scale = MatScale_Nest; 1213*a061e289SJed Brown A->ops->shift = MatShift_Nest; 1214c8883902SJed Brown 1215c8883902SJed Brown A->spptr = 0; 1216c8883902SJed Brown A->same_nonzero = PETSC_FALSE; 1217c8883902SJed Brown A->assembled = PETSC_FALSE; 1218c8883902SJed Brown 1219c8883902SJed Brown /* expose Nest api's */ 1220c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "MatNestGetSubMat_Nest", MatNestGetSubMat_Nest);CHKERRQ(ierr); 12210782ca92SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "MatNestSetSubMat_Nest", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1222c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "MatNestGetSubMats_Nest", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1223c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "MatNestGetSize_Nest", MatNestGetSize_Nest);CHKERRQ(ierr); 1224c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "MatNestSetVecType_Nest", MatNestSetVecType_Nest);CHKERRQ(ierr); 1225c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "MatNestSetSubMats_Nest", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1226c8883902SJed Brown 1227c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1228c8883902SJed Brown PetscFunctionReturn(0); 1229c8883902SJed Brown } 123086e80854SHong Zhang EXTERN_C_END 1231