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); 78336d21e7SJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&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); 198900e7ff2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetISs_C", "",0);CHKERRQ(ierr); 199900e7ff2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetLocalISs_C", "",0);CHKERRQ(ierr); 200c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "",0);CHKERRQ(ierr); 201c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "",0);CHKERRQ(ierr); 202d8588912SDave May PetscFunctionReturn(0); 203d8588912SDave May } 204d8588912SDave May 205d8588912SDave May #undef __FUNCT__ 206d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest" 207207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 208d8588912SDave May { 209d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 210d8588912SDave May PetscInt i,j; 211d8588912SDave May PetscErrorCode ierr; 212d8588912SDave May 213d8588912SDave May PetscFunctionBegin; 214d8588912SDave May for (i=0; i<vs->nr; i++) { 215d8588912SDave May for (j=0; j<vs->nc; j++) { 216d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); } 217d8588912SDave May } 218d8588912SDave May } 219d8588912SDave May PetscFunctionReturn(0); 220d8588912SDave May } 221d8588912SDave May 222d8588912SDave May #undef __FUNCT__ 223d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 224207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 225d8588912SDave May { 226d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 227d8588912SDave May PetscInt i,j; 228d8588912SDave May PetscErrorCode ierr; 229d8588912SDave May 230d8588912SDave May PetscFunctionBegin; 231d8588912SDave May for (i=0; i<vs->nr; i++) { 232d8588912SDave May for (j=0; j<vs->nc; j++) { 233d8588912SDave May if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); } 234d8588912SDave May } 235d8588912SDave May } 236d8588912SDave May PetscFunctionReturn(0); 237d8588912SDave May } 238d8588912SDave May 239d8588912SDave May #undef __FUNCT__ 240f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 241f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 242d8588912SDave May { 243207556f9SJed Brown PetscErrorCode ierr; 244f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 245f349c1fdSJed Brown PetscInt j; 246f349c1fdSJed Brown Mat sub; 247d8588912SDave May 248d8588912SDave May PetscFunctionBegin; 2498188e55aSJed Brown sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */ 250f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 2514994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 252f349c1fdSJed Brown *B = sub; 253f349c1fdSJed Brown PetscFunctionReturn(0); 254d8588912SDave May } 255d8588912SDave May 256f349c1fdSJed Brown #undef __FUNCT__ 257f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 258f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 259f349c1fdSJed Brown { 260207556f9SJed Brown PetscErrorCode ierr; 261f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 262f349c1fdSJed Brown PetscInt i; 263f349c1fdSJed Brown Mat sub; 264f349c1fdSJed Brown 265f349c1fdSJed Brown PetscFunctionBegin; 2668188e55aSJed Brown sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */ 267f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 2684994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 269f349c1fdSJed Brown *B = sub; 270f349c1fdSJed Brown PetscFunctionReturn(0); 271d8588912SDave May } 272d8588912SDave May 273f349c1fdSJed Brown #undef __FUNCT__ 274f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 275f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 276f349c1fdSJed Brown { 277f349c1fdSJed Brown PetscErrorCode ierr; 278f349c1fdSJed Brown PetscInt i; 279f349c1fdSJed Brown PetscBool flg; 280f349c1fdSJed Brown 281f349c1fdSJed Brown PetscFunctionBegin; 282f349c1fdSJed Brown PetscValidPointer(list,3); 283f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 284f349c1fdSJed Brown PetscValidIntPointer(found,5); 285f349c1fdSJed Brown *found = -1; 286f349c1fdSJed Brown for (i=0; i<n; i++) { 287207556f9SJed Brown if (!list[i]) continue; 288f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 289f349c1fdSJed Brown if (flg) { 290f349c1fdSJed Brown *found = i; 291f349c1fdSJed Brown PetscFunctionReturn(0); 292f349c1fdSJed Brown } 293f349c1fdSJed Brown } 294f349c1fdSJed Brown SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set"); 295f349c1fdSJed Brown PetscFunctionReturn(0); 296f349c1fdSJed Brown } 297f349c1fdSJed Brown 298f349c1fdSJed Brown #undef __FUNCT__ 2998188e55aSJed Brown #define __FUNCT__ "MatNestGetRow" 3008188e55aSJed Brown /* Get a block row as a new MatNest */ 3018188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 3028188e55aSJed Brown { 3038188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3048188e55aSJed Brown char keyname[256]; 3058188e55aSJed Brown PetscErrorCode ierr; 3068188e55aSJed Brown 3078188e55aSJed Brown PetscFunctionBegin; 3088188e55aSJed Brown *B = PETSC_NULL; 3098188e55aSJed Brown ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr); 3108188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 3118188e55aSJed Brown if (*B) PetscFunctionReturn(0); 3128188e55aSJed Brown 3138188e55aSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 3148188e55aSJed Brown (*B)->assembled = A->assembled; 3158188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 3168188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 3178188e55aSJed Brown PetscFunctionReturn(0); 3188188e55aSJed Brown } 3198188e55aSJed Brown 3208188e55aSJed Brown #undef __FUNCT__ 321f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 322f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 323f349c1fdSJed Brown { 324f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3258188e55aSJed Brown PetscErrorCode ierr; 3266b3a5b13SJed Brown PetscInt row,col; 3278188e55aSJed Brown PetscBool same,isFullCol; 328f349c1fdSJed Brown 329f349c1fdSJed Brown PetscFunctionBegin; 3308188e55aSJed Brown /* Check if full column space. This is a hack */ 3318188e55aSJed Brown isFullCol = PETSC_FALSE; 3328188e55aSJed Brown ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 3338188e55aSJed Brown if (same) { 33477019fcaSJed Brown PetscInt n,first,step,i,an,am,afirst,astep; 3358188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 3368188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 33777019fcaSJed Brown isFullCol = PETSC_TRUE; 33877019fcaSJed Brown for (i=0,an=0; i<vs->nc; i++) { 33977019fcaSJed Brown ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr); 34077019fcaSJed Brown ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr); 34177019fcaSJed Brown if (afirst != an || astep != step) isFullCol = PETSC_FALSE; 34277019fcaSJed Brown an += am; 34377019fcaSJed Brown } 34477019fcaSJed Brown if (an != n) isFullCol = PETSC_FALSE; 3458188e55aSJed Brown } 3468188e55aSJed Brown 3478188e55aSJed Brown if (isFullCol) { 3488188e55aSJed Brown PetscInt row; 3498188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 3508188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 3518188e55aSJed Brown } else { 352f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 353f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 354f349c1fdSJed Brown *B = vs->m[row][col]; 3558188e55aSJed Brown } 356f349c1fdSJed Brown PetscFunctionReturn(0); 357f349c1fdSJed Brown } 358f349c1fdSJed Brown 359f349c1fdSJed Brown #undef __FUNCT__ 360f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 361f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 362f349c1fdSJed Brown { 363f349c1fdSJed Brown PetscErrorCode ierr; 364f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 365f349c1fdSJed Brown Mat sub; 366f349c1fdSJed Brown 367f349c1fdSJed Brown PetscFunctionBegin; 368f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 369f349c1fdSJed Brown switch (reuse) { 370f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 3717874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 372f349c1fdSJed Brown *B = sub; 373f349c1fdSJed Brown break; 374f349c1fdSJed Brown case MAT_REUSE_MATRIX: 375f349c1fdSJed Brown if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 376f349c1fdSJed Brown break; 377f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 378f349c1fdSJed Brown break; 379f349c1fdSJed Brown } 380f349c1fdSJed Brown PetscFunctionReturn(0); 381f349c1fdSJed Brown } 382f349c1fdSJed Brown 383f349c1fdSJed Brown #undef __FUNCT__ 384f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 385f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 386f349c1fdSJed Brown { 387f349c1fdSJed Brown PetscErrorCode ierr; 388f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 389f349c1fdSJed Brown Mat sub; 390f349c1fdSJed Brown 391f349c1fdSJed Brown PetscFunctionBegin; 392f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 393f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 394f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 395f349c1fdSJed Brown *B = sub; 396d8588912SDave May PetscFunctionReturn(0); 397d8588912SDave May } 398d8588912SDave May 399d8588912SDave May #undef __FUNCT__ 400d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 401207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 402d8588912SDave May { 403d8588912SDave May PetscErrorCode ierr; 404f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 405f349c1fdSJed Brown Mat sub; 406d8588912SDave May 407d8588912SDave May PetscFunctionBegin; 408f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 409f349c1fdSJed Brown if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 410f349c1fdSJed Brown if (sub) { 411f349c1fdSJed Brown if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4126bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 413d8588912SDave May } 414d8588912SDave May PetscFunctionReturn(0); 415d8588912SDave May } 416d8588912SDave May 417d8588912SDave May #undef __FUNCT__ 4187874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 4197874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 4207874fa86SDave May { 4217874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4227874fa86SDave May PetscInt i; 4237874fa86SDave May PetscErrorCode ierr; 4247874fa86SDave May 4257874fa86SDave May PetscFunctionBegin; 4267874fa86SDave May for (i=0; i<bA->nr; i++) { 427429bac76SJed Brown Vec bv; 428429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4297874fa86SDave May if (bA->m[i][i]) { 430429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 4317874fa86SDave May } else { 432429bac76SJed Brown ierr = VecSet(bv,1.0);CHKERRQ(ierr); 4337874fa86SDave May } 434429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4357874fa86SDave May } 4367874fa86SDave May PetscFunctionReturn(0); 4377874fa86SDave May } 4387874fa86SDave May 4397874fa86SDave May #undef __FUNCT__ 4407874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 4417874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 4427874fa86SDave May { 4437874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 444429bac76SJed Brown Vec bl,*br; 4457874fa86SDave May PetscInt i,j; 4467874fa86SDave May PetscErrorCode ierr; 4477874fa86SDave May 4487874fa86SDave May PetscFunctionBegin; 449429bac76SJed Brown ierr = PetscMalloc(bA->nc*sizeof(Vec),&br);CHKERRQ(ierr); 450429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 4517874fa86SDave May for (i=0; i<bA->nr; i++) { 452429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4537874fa86SDave May for (j=0; j<bA->nc; j++) { 4547874fa86SDave May if (bA->m[i][j]) { 455429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 4567874fa86SDave May } 4577874fa86SDave May } 458a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4597874fa86SDave May } 460429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 461429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 4627874fa86SDave May PetscFunctionReturn(0); 4637874fa86SDave May } 4647874fa86SDave May 4657874fa86SDave May #undef __FUNCT__ 466a061e289SJed Brown #define __FUNCT__ "MatScale_Nest" 467a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 468a061e289SJed Brown { 469a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 470a061e289SJed Brown PetscInt i,j; 471a061e289SJed Brown PetscErrorCode ierr; 472a061e289SJed Brown 473a061e289SJed Brown PetscFunctionBegin; 474a061e289SJed Brown for (i=0; i<bA->nr; i++) { 475a061e289SJed Brown for (j=0; j<bA->nc; j++) { 476a061e289SJed Brown if (bA->m[i][j]) { 477a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 478a061e289SJed Brown } 479a061e289SJed Brown } 480a061e289SJed Brown } 481a061e289SJed Brown PetscFunctionReturn(0); 482a061e289SJed Brown } 483a061e289SJed Brown 484a061e289SJed Brown #undef __FUNCT__ 485a061e289SJed Brown #define __FUNCT__ "MatShift_Nest" 486a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 487a061e289SJed Brown { 488a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 489a061e289SJed Brown PetscInt i; 490a061e289SJed Brown PetscErrorCode ierr; 491a061e289SJed Brown 492a061e289SJed Brown PetscFunctionBegin; 493a061e289SJed Brown for (i=0; i<bA->nr; i++) { 494a061e289SJed 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); 495a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 496a061e289SJed Brown } 497a061e289SJed Brown PetscFunctionReturn(0); 498a061e289SJed Brown } 499a061e289SJed Brown 500a061e289SJed Brown #undef __FUNCT__ 501d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest" 502207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 503d8588912SDave May { 504d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 505d8588912SDave May Vec *L,*R; 506d8588912SDave May MPI_Comm comm; 507d8588912SDave May PetscInt i,j; 508d8588912SDave May PetscErrorCode ierr; 509d8588912SDave May 510d8588912SDave May PetscFunctionBegin; 511d8588912SDave May comm = ((PetscObject)A)->comm; 512d8588912SDave May if (right) { 513d8588912SDave May /* allocate R */ 514d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr); 515d8588912SDave May /* Create the right vectors */ 516d8588912SDave May for (j=0; j<bA->nc; j++) { 517d8588912SDave May for (i=0; i<bA->nr; i++) { 518d8588912SDave May if (bA->m[i][j]) { 519d8588912SDave May ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 520d8588912SDave May break; 521d8588912SDave May } 522d8588912SDave May } 523d8588912SDave May if (i==bA->nr) { 524d8588912SDave May /* have an empty column */ 525d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 526d8588912SDave May } 527d8588912SDave May } 528f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 529d8588912SDave May /* hand back control to the nest vector */ 530d8588912SDave May for (j=0; j<bA->nc; j++) { 5316bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 532d8588912SDave May } 533d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 534d8588912SDave May } 535d8588912SDave May 536d8588912SDave May if (left) { 537d8588912SDave May /* allocate L */ 538d8588912SDave May ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr); 539d8588912SDave May /* Create the left vectors */ 540d8588912SDave May for (i=0; i<bA->nr; i++) { 541d8588912SDave May for (j=0; j<bA->nc; j++) { 542d8588912SDave May if (bA->m[i][j]) { 543d8588912SDave May ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 544d8588912SDave May break; 545d8588912SDave May } 546d8588912SDave May } 547d8588912SDave May if (j==bA->nc) { 548d8588912SDave May /* have an empty row */ 549d8588912SDave May SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 550d8588912SDave May } 551d8588912SDave May } 552d8588912SDave May 553f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 554d8588912SDave May for (i=0; i<bA->nr; i++) { 5556bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 556d8588912SDave May } 557d8588912SDave May 558d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 559d8588912SDave May } 560d8588912SDave May PetscFunctionReturn(0); 561d8588912SDave May } 562d8588912SDave May 563d8588912SDave May #undef __FUNCT__ 564d8588912SDave May #define __FUNCT__ "MatView_Nest" 565207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 566d8588912SDave May { 567d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 568d8588912SDave May PetscBool isascii; 569d8588912SDave May PetscInt i,j; 570d8588912SDave May PetscErrorCode ierr; 571d8588912SDave May 572d8588912SDave May PetscFunctionBegin; 573d8588912SDave May ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 574d8588912SDave May if (isascii) { 575d8588912SDave May 576d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n" ); 577d8588912SDave May PetscViewerASCIIPushTab( viewer ); /* push0 */ 578d8588912SDave May PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 579d8588912SDave May 580d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" ); 581d8588912SDave May for (i=0; i<bA->nr; i++) { 582d8588912SDave May for (j=0; j<bA->nc; j++) { 583d8588912SDave May const MatType type; 584270f95d7SJed Brown char name[256] = "",prefix[256] = ""; 585d8588912SDave May PetscInt NR,NC; 586d8588912SDave May PetscBool isNest = PETSC_FALSE; 587d8588912SDave May 588d8588912SDave May if (!bA->m[i][j]) { 589d8588912SDave May PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j); 590d8588912SDave May continue; 591d8588912SDave May } 592d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 593d8588912SDave May ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr); 594270f95d7SJed Brown if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof name,"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 595270f95d7SJed Brown if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof prefix,"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 596d8588912SDave May ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 597d8588912SDave May 598270f95d7SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 599d8588912SDave May 600d8588912SDave May if (isNest) { 601270f95d7SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 602d8588912SDave May ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 603270f95d7SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 604d8588912SDave May } 605d8588912SDave May } 606d8588912SDave May } 607d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 608d8588912SDave May } 609d8588912SDave May PetscFunctionReturn(0); 610d8588912SDave May } 611d8588912SDave May 612d8588912SDave May #undef __FUNCT__ 613d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 614207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 615d8588912SDave May { 616d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 617d8588912SDave May PetscInt i,j; 618d8588912SDave May PetscErrorCode ierr; 619d8588912SDave May 620d8588912SDave May PetscFunctionBegin; 621d8588912SDave May for (i=0; i<bA->nr; i++) { 622d8588912SDave May for (j=0; j<bA->nc; j++) { 623d8588912SDave May if (!bA->m[i][j]) continue; 624d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 625d8588912SDave May } 626d8588912SDave May } 627d8588912SDave May PetscFunctionReturn(0); 628d8588912SDave May } 629d8588912SDave May 630d8588912SDave May #undef __FUNCT__ 631d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 632207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 633d8588912SDave May { 634d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 635841e96a3SJed Brown Mat *b; 636841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 637d8588912SDave May PetscErrorCode ierr; 638d8588912SDave May 639d8588912SDave May PetscFunctionBegin; 640841e96a3SJed Brown ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 641841e96a3SJed Brown for (i=0; i<nr; i++) { 642841e96a3SJed Brown for (j=0; j<nc; j++) { 643841e96a3SJed Brown if (bA->m[i][j]) { 644841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 645841e96a3SJed Brown } else { 646841e96a3SJed Brown b[i*nc+j] = PETSC_NULL; 647d8588912SDave May } 648d8588912SDave May } 649d8588912SDave May } 650f349c1fdSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 651841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 652841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 6536bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 654d8588912SDave May } 655d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 656d8588912SDave May 657841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 658841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 659d8588912SDave May PetscFunctionReturn(0); 660d8588912SDave May } 661d8588912SDave May 662d8588912SDave May /* nest api */ 663d8588912SDave May EXTERN_C_BEGIN 664d8588912SDave May #undef __FUNCT__ 665d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 666d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 667d8588912SDave May { 668d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 669d8588912SDave May PetscFunctionBegin; 670d8588912SDave May if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 671d8588912SDave May if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 672d8588912SDave May *mat = bA->m[idxm][jdxm]; 673d8588912SDave May PetscFunctionReturn(0); 674d8588912SDave May } 675d8588912SDave May EXTERN_C_END 676d8588912SDave May 677d8588912SDave May #undef __FUNCT__ 678d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 679d8588912SDave May /*@C 680d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 681d8588912SDave May 682d8588912SDave May Not collective 683d8588912SDave May 684d8588912SDave May Input Parameters: 685629881c0SJed Brown + A - nest matrix 686d8588912SDave May . idxm - index of the matrix within the nest matrix 687629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 688d8588912SDave May 689d8588912SDave May Output Parameter: 690d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 691d8588912SDave May 692d8588912SDave May Level: developer 693d8588912SDave May 694d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 695d8588912SDave May @*/ 6967087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 697d8588912SDave May { 698699a902aSJed Brown PetscErrorCode ierr; 699d8588912SDave May 700d8588912SDave May PetscFunctionBegin; 701699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 702d8588912SDave May PetscFunctionReturn(0); 703d8588912SDave May } 704d8588912SDave May 705d8588912SDave May EXTERN_C_BEGIN 706d8588912SDave May #undef __FUNCT__ 7070782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 7080782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 7090782ca92SJed Brown { 7100782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 7110782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 7120782ca92SJed Brown PetscErrorCode ierr; 7130782ca92SJed Brown 7140782ca92SJed Brown PetscFunctionBegin; 7150782ca92SJed Brown if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 7160782ca92SJed Brown if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 7170782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 7180782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 7190782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 7200782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 7210782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 7220782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 7230782ca92SJed 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); 7240782ca92SJed 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); 7250782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 7260782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 7270782ca92SJed Brown bA->m[idxm][jdxm] = mat; 7280782ca92SJed Brown PetscFunctionReturn(0); 7290782ca92SJed Brown } 7300782ca92SJed Brown EXTERN_C_END 7310782ca92SJed Brown 7320782ca92SJed Brown #undef __FUNCT__ 7330782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 7340782ca92SJed Brown /*@C 7350782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 7360782ca92SJed Brown 7370782ca92SJed Brown Logically collective on the submatrix communicator 7380782ca92SJed Brown 7390782ca92SJed Brown Input Parameters: 7400782ca92SJed Brown + A - nest matrix 7410782ca92SJed Brown . idxm - index of the matrix within the nest matrix 7420782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 7430782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 7440782ca92SJed Brown 7450782ca92SJed Brown Notes: 7460782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 7470782ca92SJed Brown 7480782ca92SJed Brown This increments the reference count of the submatrix. 7490782ca92SJed Brown 7500782ca92SJed Brown Level: developer 7510782ca92SJed Brown 7520782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 7530782ca92SJed Brown @*/ 7540782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 7550782ca92SJed Brown { 7560782ca92SJed Brown PetscErrorCode ierr; 7570782ca92SJed Brown 7580782ca92SJed Brown PetscFunctionBegin; 7590782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 7600782ca92SJed Brown PetscFunctionReturn(0); 7610782ca92SJed Brown } 7620782ca92SJed Brown 7630782ca92SJed Brown EXTERN_C_BEGIN 7640782ca92SJed Brown #undef __FUNCT__ 765d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 766d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 767d8588912SDave May { 768d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 769d8588912SDave May PetscFunctionBegin; 770d8588912SDave May if (M) { *M = bA->nr; } 771d8588912SDave May if (N) { *N = bA->nc; } 772d8588912SDave May if (mat) { *mat = bA->m; } 773d8588912SDave May PetscFunctionReturn(0); 774d8588912SDave May } 775d8588912SDave May EXTERN_C_END 776d8588912SDave May 777d8588912SDave May #undef __FUNCT__ 778d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 779d8588912SDave May /*@C 780d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 781d8588912SDave May 782d8588912SDave May Not collective 783d8588912SDave May 784d8588912SDave May Input Parameters: 785629881c0SJed Brown . A - nest matrix 786d8588912SDave May 787d8588912SDave May Output Parameter: 788629881c0SJed Brown + M - number of rows in the nest matrix 789d8588912SDave May . N - number of cols in the nest matrix 790629881c0SJed Brown - mat - 2d array of matrices 791d8588912SDave May 792d8588912SDave May Notes: 793d8588912SDave May 794d8588912SDave May The user should not free the array mat. 795d8588912SDave May 796d8588912SDave May Level: developer 797d8588912SDave May 798d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 799d8588912SDave May @*/ 8007087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 801d8588912SDave May { 802699a902aSJed Brown PetscErrorCode ierr; 803d8588912SDave May 804d8588912SDave May PetscFunctionBegin; 805699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 806d8588912SDave May PetscFunctionReturn(0); 807d8588912SDave May } 808d8588912SDave May 809d8588912SDave May EXTERN_C_BEGIN 810d8588912SDave May #undef __FUNCT__ 811d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 8127087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 813d8588912SDave May { 814d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 815d8588912SDave May 816d8588912SDave May PetscFunctionBegin; 817d8588912SDave May if (M) { *M = bA->nr; } 818d8588912SDave May if (N) { *N = bA->nc; } 819d8588912SDave May PetscFunctionReturn(0); 820d8588912SDave May } 821d8588912SDave May EXTERN_C_END 822d8588912SDave May 823d8588912SDave May #undef __FUNCT__ 824d8588912SDave May #define __FUNCT__ "MatNestGetSize" 825d8588912SDave May /*@C 826d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 827d8588912SDave May 828d8588912SDave May Not collective 829d8588912SDave May 830d8588912SDave May Input Parameters: 831d8588912SDave May . A - nest matrix 832d8588912SDave May 833d8588912SDave May Output Parameter: 834629881c0SJed Brown + M - number of rows in the nested mat 835629881c0SJed Brown - N - number of cols in the nested mat 836d8588912SDave May 837d8588912SDave May Notes: 838d8588912SDave May 839d8588912SDave May Level: developer 840d8588912SDave May 841d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 842d8588912SDave May @*/ 8437087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 844d8588912SDave May { 845699a902aSJed Brown PetscErrorCode ierr; 846d8588912SDave May 847d8588912SDave May PetscFunctionBegin; 848699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 849d8588912SDave May PetscFunctionReturn(0); 850d8588912SDave May } 851d8588912SDave May 852900e7ff2SJed Brown #undef __FUNCT__ 853900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest" 854900e7ff2SJed Brown PETSC_EXTERN_C PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 855900e7ff2SJed Brown { 856900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 857900e7ff2SJed Brown PetscInt i; 858900e7ff2SJed Brown 859900e7ff2SJed Brown PetscFunctionBegin; 860900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 861900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 862900e7ff2SJed Brown PetscFunctionReturn(0); 863900e7ff2SJed Brown } 864900e7ff2SJed Brown 865900e7ff2SJed Brown #undef __FUNCT__ 866900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs" 867900e7ff2SJed Brown /*@C 868900e7ff2SJed Brown MatNestGetISs - Returns the index sets partitioning the row and column spaces 869900e7ff2SJed Brown 870900e7ff2SJed Brown Not collective 871900e7ff2SJed Brown 872900e7ff2SJed Brown Input Parameters: 873900e7ff2SJed Brown . A - nest matrix 874900e7ff2SJed Brown 875900e7ff2SJed Brown Output Parameter: 876900e7ff2SJed Brown + rows - array of row index sets 877900e7ff2SJed Brown - cols - array of column index sets 878900e7ff2SJed Brown 879900e7ff2SJed Brown Level: advanced 880900e7ff2SJed Brown 881900e7ff2SJed Brown Notes: 882900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 883900e7ff2SJed Brown 884900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs() 885900e7ff2SJed Brown @*/ 886900e7ff2SJed Brown PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 887900e7ff2SJed Brown { 888900e7ff2SJed Brown PetscErrorCode ierr; 889900e7ff2SJed Brown 890900e7ff2SJed Brown PetscFunctionBegin; 891900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 892900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 893900e7ff2SJed Brown PetscFunctionReturn(0); 894900e7ff2SJed Brown } 895900e7ff2SJed Brown 896900e7ff2SJed Brown #undef __FUNCT__ 897900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest" 898900e7ff2SJed Brown PETSC_EXTERN_C PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 899900e7ff2SJed Brown { 900900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 901900e7ff2SJed Brown PetscInt i; 902900e7ff2SJed Brown 903900e7ff2SJed Brown PetscFunctionBegin; 904900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 905900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 906900e7ff2SJed Brown PetscFunctionReturn(0); 907900e7ff2SJed Brown } 908900e7ff2SJed Brown 909900e7ff2SJed Brown #undef __FUNCT__ 910900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs" 911900e7ff2SJed Brown /*@C 912900e7ff2SJed Brown MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 913900e7ff2SJed Brown 914900e7ff2SJed Brown Not collective 915900e7ff2SJed Brown 916900e7ff2SJed Brown Input Parameters: 917900e7ff2SJed Brown . A - nest matrix 918900e7ff2SJed Brown 919900e7ff2SJed Brown Output Parameter: 920900e7ff2SJed Brown + rows - array of row index sets (or PETSC_NULL to ignore) 921900e7ff2SJed Brown - cols - array of column index sets (or PETSC_NULL to ignore) 922900e7ff2SJed Brown 923900e7ff2SJed Brown Level: advanced 924900e7ff2SJed Brown 925900e7ff2SJed Brown Notes: 926900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 927900e7ff2SJed Brown 928900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs() 929900e7ff2SJed Brown @*/ 930900e7ff2SJed Brown PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 931900e7ff2SJed Brown { 932900e7ff2SJed Brown PetscErrorCode ierr; 933900e7ff2SJed Brown 934900e7ff2SJed Brown PetscFunctionBegin; 935900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 936900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 937900e7ff2SJed Brown PetscFunctionReturn(0); 938900e7ff2SJed Brown } 939900e7ff2SJed Brown 940207556f9SJed Brown EXTERN_C_BEGIN 941207556f9SJed Brown #undef __FUNCT__ 942207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 9437087cfbeSBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,const VecType vtype) 944207556f9SJed Brown { 945207556f9SJed Brown PetscErrorCode ierr; 946207556f9SJed Brown PetscBool flg; 947207556f9SJed Brown 948207556f9SJed Brown PetscFunctionBegin; 949207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 950207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 951207556f9SJed Brown A->ops->getvecs = flg ? MatGetVecs_Nest : 0; 952207556f9SJed Brown PetscFunctionReturn(0); 953207556f9SJed Brown } 954207556f9SJed Brown EXTERN_C_END 955207556f9SJed Brown 956207556f9SJed Brown #undef __FUNCT__ 957207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 958207556f9SJed Brown /*@C 959207556f9SJed Brown MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 960207556f9SJed Brown 961207556f9SJed Brown Not collective 962207556f9SJed Brown 963207556f9SJed Brown Input Parameters: 964207556f9SJed Brown + A - nest matrix 965207556f9SJed Brown - vtype - type to use for creating vectors 966207556f9SJed Brown 967207556f9SJed Brown Notes: 968207556f9SJed Brown 969207556f9SJed Brown Level: developer 970207556f9SJed Brown 971207556f9SJed Brown .seealso: MatGetVecs() 972207556f9SJed Brown @*/ 9737087cfbeSBarry Smith PetscErrorCode MatNestSetVecType(Mat A,const VecType vtype) 974207556f9SJed Brown { 975207556f9SJed Brown PetscErrorCode ierr; 976207556f9SJed Brown 977207556f9SJed Brown PetscFunctionBegin; 978207556f9SJed Brown ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr); 979207556f9SJed Brown PetscFunctionReturn(0); 980207556f9SJed Brown } 981207556f9SJed Brown 982c8883902SJed Brown EXTERN_C_BEGIN 983d8588912SDave May #undef __FUNCT__ 984c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 985c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 986d8588912SDave May { 987c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 988c8883902SJed Brown PetscInt i,j,m,n,M,N; 989d8588912SDave May PetscErrorCode ierr; 990d8588912SDave May 991d8588912SDave May PetscFunctionBegin; 992c8883902SJed Brown s->nr = nr; 993c8883902SJed Brown s->nc = nc; 994d8588912SDave May 995c8883902SJed Brown /* Create space for submatrices */ 996c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 997c8883902SJed Brown for (i=0; i<nr; i++) { 998c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 999d8588912SDave May } 1000c8883902SJed Brown for (i=0; i<nr; i++) { 1001c8883902SJed Brown for (j=0; j<nc; j++) { 1002c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1003c8883902SJed Brown if (a[i*nc+j]) { 1004c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1005d8588912SDave May } 1006d8588912SDave May } 1007d8588912SDave May } 1008d8588912SDave May 10098188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1010d8588912SDave May 1011c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 1012c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 1013c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1014c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1015d8588912SDave May 10168188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1017d8588912SDave May 1018c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1019c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1020c8883902SJed Brown ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 1021c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1022c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1023c8883902SJed Brown ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 1024c8883902SJed Brown 1025c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1026c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1027c8883902SJed Brown 1028c8883902SJed Brown ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 1029c8883902SJed Brown ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr); 1030c8883902SJed Brown ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr); 1031d8588912SDave May PetscFunctionReturn(0); 1032d8588912SDave May } 1033c8883902SJed Brown EXTERN_C_END 1034d8588912SDave May 1035c8883902SJed Brown #undef __FUNCT__ 1036c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 1037c8883902SJed Brown /*@ 1038c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1039c8883902SJed Brown 1040c8883902SJed Brown Collective on Mat 1041c8883902SJed Brown 1042c8883902SJed Brown Input Parameter: 1043c8883902SJed Brown + N - nested matrix 1044c8883902SJed Brown . nr - number of nested row blocks 1045c8883902SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1046c8883902SJed Brown . nc - number of nested column blocks 1047c8883902SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1048c8883902SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1049c8883902SJed Brown 1050c8883902SJed Brown Level: advanced 1051c8883902SJed Brown 1052c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 1053c8883902SJed Brown @*/ 1054c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1055c8883902SJed Brown { 1056c8883902SJed Brown PetscErrorCode ierr; 1057c8883902SJed Brown PetscInt i; 1058c8883902SJed Brown 1059c8883902SJed Brown PetscFunctionBegin; 1060c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1061c8883902SJed Brown if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1062c8883902SJed Brown if (nr && is_row) { 1063c8883902SJed Brown PetscValidPointer(is_row,3); 1064c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1065c8883902SJed Brown } 1066c8883902SJed Brown if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 10671664e352SJed Brown if (nc && is_col) { 1068c8883902SJed Brown PetscValidPointer(is_col,5); 1069c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1070c8883902SJed Brown } 1071c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 1072c8883902SJed 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); 1073c8883902SJed Brown PetscFunctionReturn(0); 1074c8883902SJed Brown } 1075d8588912SDave May 107677019fcaSJed Brown #undef __FUNCT__ 107777019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private" 107877019fcaSJed Brown static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb) 107977019fcaSJed Brown { 108077019fcaSJed Brown PetscErrorCode ierr; 108177019fcaSJed Brown PetscBool flg; 108277019fcaSJed Brown PetscInt i,j,m,mi,*ix; 108377019fcaSJed Brown 108477019fcaSJed Brown PetscFunctionBegin; 108577019fcaSJed Brown for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 108677019fcaSJed Brown if (islocal[i]) { 108777019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 108877019fcaSJed Brown flg = PETSC_TRUE; /* We found a non-trivial entry */ 108977019fcaSJed Brown } else { 109077019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 109177019fcaSJed Brown } 109277019fcaSJed Brown m += mi; 109377019fcaSJed Brown } 109477019fcaSJed Brown if (flg) { 109577019fcaSJed Brown ierr = PetscMalloc(m*sizeof(*ix),&ix);CHKERRQ(ierr); 109677019fcaSJed Brown for (i=0,n=0; i<n; i++) { 109777019fcaSJed Brown ISLocalToGlobalMapping smap = PETSC_NULL; 109877019fcaSJed Brown VecScatter scat; 109977019fcaSJed Brown IS isreq; 110077019fcaSJed Brown Vec lvec,gvec; 11013361c9a7SJed Brown union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x; 110277019fcaSJed Brown Mat sub; 110377019fcaSJed Brown 110477019fcaSJed Brown if (sizeof (*x) != sizeof(PetscScalar)) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support when scalars smaller than integers"); 110577019fcaSJed Brown if (colflg) { 110677019fcaSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 110777019fcaSJed Brown } else { 110877019fcaSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 110977019fcaSJed Brown } 111077019fcaSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,PETSC_NULL);CHKERRQ(ierr);} 111177019fcaSJed Brown if (islocal[i]) { 111277019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 111377019fcaSJed Brown } else { 111477019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 111577019fcaSJed Brown } 111677019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = j; 111777019fcaSJed Brown if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);} 111877019fcaSJed Brown /* 111977019fcaSJed Brown Now we need to extract the monolithic global indices that correspond to the given split global indices. 112077019fcaSJed 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. 112177019fcaSJed Brown The approach here is ugly because it uses VecScatter to move indices. 112277019fcaSJed Brown */ 112377019fcaSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr); 112477019fcaSJed Brown ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr); 112577019fcaSJed Brown ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr); 112677019fcaSJed Brown ierr = VecScatterCreate(gvec,isreq,lvec,PETSC_NULL,&scat);CHKERRQ(ierr); 112777019fcaSJed Brown ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 112877019fcaSJed Brown for (j=0; j<mi; j++) x[j].integer = ix[m+j]; 112977019fcaSJed Brown ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 113077019fcaSJed Brown ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 113177019fcaSJed Brown ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 113277019fcaSJed Brown ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 113377019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = x[j].integer; 113477019fcaSJed Brown ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 113577019fcaSJed Brown ierr = VecDestroy(&lvec);CHKERRQ(ierr); 113677019fcaSJed Brown ierr = VecDestroy(&gvec);CHKERRQ(ierr); 113777019fcaSJed Brown ierr = ISDestroy(&isreq);CHKERRQ(ierr); 113877019fcaSJed Brown ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 113977019fcaSJed Brown m += mi; 114077019fcaSJed Brown } 114177019fcaSJed Brown ierr = ISLocalToGlobalMappingCreate(((PetscObject)A)->comm,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 114277019fcaSJed Brown *ltogb = PETSC_NULL; 114377019fcaSJed Brown } else { 114477019fcaSJed Brown *ltog = PETSC_NULL; 114577019fcaSJed Brown *ltogb = PETSC_NULL; 114677019fcaSJed Brown } 114777019fcaSJed Brown PetscFunctionReturn(0); 114877019fcaSJed Brown } 114977019fcaSJed Brown 115077019fcaSJed Brown 1151d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1152d8588912SDave May /* 1153d8588912SDave May nprocessors = NP 1154d8588912SDave May Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) ) 1155d8588912SDave May proc 0: => (g_0,h_0,) 1156d8588912SDave May proc 1: => (g_1,h_1,) 1157d8588912SDave May ... 1158d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1159d8588912SDave May 1160d8588912SDave May proc 0: proc 1: proc nprocs-1: 1161d8588912SDave May is[0] = ( 0,1,2,...,nlocal(g_0)-1 ) ( 0,1,...,nlocal(g_1)-1 ) ( 0,1,...,nlocal(g_NP-1) ) 1162d8588912SDave May 1163d8588912SDave May proc 0: 1164d8588912SDave May is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 ) 1165d8588912SDave May proc 1: 1166d8588912SDave May is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 ) 1167d8588912SDave May 1168d8588912SDave May proc NP-1: 1169d8588912SDave May is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 ) 1170d8588912SDave May */ 1171d8588912SDave May #undef __FUNCT__ 1172d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1173841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1174d8588912SDave May { 1175e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 11768188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1177d8588912SDave May PetscErrorCode ierr; 11782ae74bdbSJed Brown Mat sub; 1179d8588912SDave May 1180d8588912SDave May PetscFunctionBegin; 11818188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 11828188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1183d8588912SDave May if (is_row) { /* valid IS is passed in */ 1184d8588912SDave May /* refs on is[] are incremeneted */ 1185e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1186d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1187e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1188d8588912SDave May } 11892ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 11908188e55aSJed Brown nsum = 0; 11918188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 11928188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 11938188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 11948188e55aSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 11958188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 11968188e55aSJed Brown nsum += n; 11978188e55aSJed Brown } 11988188e55aSJed Brown offset = 0; 11993fda65a0SJose Roman #if defined(PETSC_HAVE_MPI_EXSCAN) 12008188e55aSJed Brown ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 12013fda65a0SJose Roman #else 12023fda65a0SJose Roman SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sorry but this code requires MPI_EXSCAN that doesn't exist on your machine's version of MPI, install a MPI2 with PETSc to get this functionality"); 12033fda65a0SJose Roman #endif 1204e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1205f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12062ae74bdbSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 12072ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1208e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1209e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 12102ae74bdbSJed Brown offset += n; 1211d8588912SDave May } 1212d8588912SDave May } 1213d8588912SDave May 1214d8588912SDave May if (is_col) { /* valid IS is passed in */ 1215d8588912SDave May /* refs on is[] are incremeneted */ 1216e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1217d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1218e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1219d8588912SDave May } 12202ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 12212ae74bdbSJed Brown offset = A->cmap->rstart; 12228188e55aSJed Brown nsum = 0; 12238188e55aSJed Brown for (j=0; j<vs->nc; j++) { 12248188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 12258188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 12268188e55aSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 12278188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12288188e55aSJed Brown nsum += n; 12298188e55aSJed Brown } 12308188e55aSJed Brown offset = 0; 1231d9a7f1c0SMark F. Adams #if defined(PETSC_HAVE_MPI_EXSCAN) 12328188e55aSJed Brown ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 1233d9a7f1c0SMark F. Adams #else 1234d9a7f1c0SMark F. Adams SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sorry but this code requires MPI_EXSCAN that doesn't exist on your machine's version of MPI, install a MPI2 with PETSc to get this functionality"); 1235d9a7f1c0SMark F. Adams #endif 1236e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1237f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 12382ae74bdbSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 12392ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1240e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1241e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 12422ae74bdbSJed Brown offset += n; 1243d8588912SDave May } 1244d8588912SDave May } 1245e2d7f03fSJed Brown 1246e2d7f03fSJed Brown /* Set up the local ISs */ 1247e2d7f03fSJed Brown ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1248e2d7f03fSJed Brown ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1249e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1250e2d7f03fSJed Brown IS isloc; 12518188e55aSJed Brown ISLocalToGlobalMapping rmap = PETSC_NULL; 1252e2d7f03fSJed Brown PetscInt nlocal,bs; 1253e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12548188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);} 1255207556f9SJed Brown if (rmap) { 1256e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1257e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1258e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1259e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1260207556f9SJed Brown } else { 1261207556f9SJed Brown nlocal = 0; 1262207556f9SJed Brown isloc = PETSC_NULL; 1263207556f9SJed Brown } 1264e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1265e2d7f03fSJed Brown offset += nlocal; 1266e2d7f03fSJed Brown } 12678188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1268e2d7f03fSJed Brown IS isloc; 12698188e55aSJed Brown ISLocalToGlobalMapping cmap = PETSC_NULL; 1270e2d7f03fSJed Brown PetscInt nlocal,bs; 1271e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 12728188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);} 1273207556f9SJed Brown if (cmap) { 1274e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1275e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1276e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1277e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1278207556f9SJed Brown } else { 1279207556f9SJed Brown nlocal = 0; 1280207556f9SJed Brown isloc = PETSC_NULL; 1281207556f9SJed Brown } 1282e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1283e2d7f03fSJed Brown offset += nlocal; 1284e2d7f03fSJed Brown } 12850189643fSJed Brown 128677019fcaSJed Brown /* Set up the aggregate ISLocalToGlobalMapping */ 128777019fcaSJed Brown { 128877019fcaSJed Brown ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb; 128977019fcaSJed Brown ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr); 129077019fcaSJed Brown ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr); 129177019fcaSJed Brown if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 129277019fcaSJed Brown if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);} 129377019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 129477019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr); 129577019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 129677019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr); 129777019fcaSJed Brown } 129877019fcaSJed Brown 12990189643fSJed Brown #if defined(PETSC_USE_DEBUG) 13000189643fSJed Brown for (i=0; i<vs->nr; i++) { 13010189643fSJed Brown for (j=0; j<vs->nc; j++) { 13020189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 13030189643fSJed Brown Mat B = vs->m[i][j]; 13040189643fSJed Brown if (!B) continue; 13050189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 13060189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 13070189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 13080189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 13090189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 13100189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 13110189643fSJed 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); 13120189643fSJed 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); 13130189643fSJed Brown } 13140189643fSJed Brown } 13150189643fSJed Brown #endif 1316a061e289SJed Brown 1317a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1318a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1319a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1320a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1321a061e289SJed Brown } 1322a061e289SJed Brown } 1323a061e289SJed Brown A->assembled = PETSC_TRUE; 1324d8588912SDave May PetscFunctionReturn(0); 1325d8588912SDave May } 1326d8588912SDave May 1327d8588912SDave May #undef __FUNCT__ 1328d8588912SDave May #define __FUNCT__ "MatCreateNest" 1329659c6bb0SJed Brown /*@ 1330659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1331659c6bb0SJed Brown 1332659c6bb0SJed Brown Collective on Mat 1333659c6bb0SJed Brown 1334659c6bb0SJed Brown Input Parameter: 1335659c6bb0SJed Brown + comm - Communicator for the new Mat 1336659c6bb0SJed Brown . nr - number of nested row blocks 1337659c6bb0SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1338659c6bb0SJed Brown . nc - number of nested column blocks 1339659c6bb0SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1340659c6bb0SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1341659c6bb0SJed Brown 1342659c6bb0SJed Brown Output Parameter: 1343659c6bb0SJed Brown . B - new matrix 1344659c6bb0SJed Brown 1345659c6bb0SJed Brown Level: advanced 1346659c6bb0SJed Brown 1347950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST 1348659c6bb0SJed Brown @*/ 13497087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1350d8588912SDave May { 1351d8588912SDave May Mat A; 1352d8588912SDave May PetscErrorCode ierr; 1353d8588912SDave May 1354d8588912SDave May PetscFunctionBegin; 1355c8883902SJed Brown *B = 0; 1356d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1357c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 1358*7ae8954aSSatish Balay ierr = MatSetUp(A);CHKERRQ(ierr); 1359c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1360d8588912SDave May *B = A; 1361d8588912SDave May PetscFunctionReturn(0); 1362d8588912SDave May } 1363659c6bb0SJed Brown 1364659c6bb0SJed Brown /*MC 1365659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1366659c6bb0SJed Brown 1367659c6bb0SJed Brown Level: intermediate 1368659c6bb0SJed Brown 1369659c6bb0SJed Brown Notes: 1370659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1371659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1372950540a4SJed Brown It is usually used with DMComposite and DMCreateMatrix() 1373659c6bb0SJed Brown 1374659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1375659c6bb0SJed Brown M*/ 1376c8883902SJed Brown EXTERN_C_BEGIN 1377c8883902SJed Brown #undef __FUNCT__ 1378c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 1379c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A) 1380c8883902SJed Brown { 1381c8883902SJed Brown Mat_Nest *s; 1382c8883902SJed Brown PetscErrorCode ierr; 1383c8883902SJed Brown 1384c8883902SJed Brown PetscFunctionBegin; 1385c8883902SJed Brown ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr); 1386c8883902SJed Brown A->data = (void*)s; 1387c8883902SJed Brown s->nr = s->nc = -1; 1388c8883902SJed Brown s->m = PETSC_NULL; 1389c8883902SJed Brown 1390c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 1391c8883902SJed Brown A->ops->mult = MatMult_Nest; 13929194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 1393c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 13949194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1395c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1396c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1397c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1398c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1399c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1400c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1401c8883902SJed Brown A->ops->view = MatView_Nest; 1402c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1403c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1404c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1405429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 1406429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 1407a061e289SJed Brown A->ops->scale = MatScale_Nest; 1408a061e289SJed Brown A->ops->shift = MatShift_Nest; 1409c8883902SJed Brown 1410c8883902SJed Brown A->spptr = 0; 1411c8883902SJed Brown A->same_nonzero = PETSC_FALSE; 1412c8883902SJed Brown A->assembled = PETSC_FALSE; 1413c8883902SJed Brown 1414c8883902SJed Brown /* expose Nest api's */ 1415c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "MatNestGetSubMat_Nest", MatNestGetSubMat_Nest);CHKERRQ(ierr); 14160782ca92SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "MatNestSetSubMat_Nest", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1417c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "MatNestGetSubMats_Nest", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1418c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "MatNestGetSize_Nest", MatNestGetSize_Nest);CHKERRQ(ierr); 1419900e7ff2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetISs_C", "MatNestGetISs_Nest", MatNestGetISs_Nest);CHKERRQ(ierr); 1420900e7ff2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetLocalISs_C", "MatNestGetLocalISs_Nest", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 1421c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "MatNestSetVecType_Nest", MatNestSetVecType_Nest);CHKERRQ(ierr); 1422c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "MatNestSetSubMats_Nest", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1423c8883902SJed Brown 1424c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1425c8883902SJed Brown PetscFunctionReturn(0); 1426c8883902SJed Brown } 142786e80854SHong Zhang EXTERN_C_END 1428