1d8588912SDave May 2ccd284c7SBarry Smith #include "../src/mat/impls/nest/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++) { 1066c75ac25SJed Brown if (!bA->m[i][j]) 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++) { 1366c75ac25SJed Brown if (!bA->m[i][j]) 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++) { 216e7c19651SJed Brown if (vs->m[i][j]) { 217e7c19651SJed Brown ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); 218e7c19651SJed Brown if (!vs->splitassembly) { 219e7c19651SJed Brown /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested 220e7c19651SJed Brown * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was 221e7c19651SJed Brown * already performing an assembly, but the result would by more complicated and appears to offer less 222e7c19651SJed Brown * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an 223e7c19651SJed Brown * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives". 224e7c19651SJed Brown */ 225e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 226e7c19651SJed Brown } 227e7c19651SJed Brown } 228d8588912SDave May } 229d8588912SDave May } 230d8588912SDave May PetscFunctionReturn(0); 231d8588912SDave May } 232d8588912SDave May 233d8588912SDave May #undef __FUNCT__ 234d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 235207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 236d8588912SDave May { 237d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 238d8588912SDave May PetscInt i,j; 239d8588912SDave May PetscErrorCode ierr; 240d8588912SDave May 241d8588912SDave May PetscFunctionBegin; 242d8588912SDave May for (i=0; i<vs->nr; i++) { 243d8588912SDave May for (j=0; j<vs->nc; j++) { 244e7c19651SJed Brown if (vs->m[i][j]) { 245e7c19651SJed Brown if (vs->splitassembly) { 246e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 247e7c19651SJed Brown } 248e7c19651SJed Brown } 249d8588912SDave May } 250d8588912SDave May } 251d8588912SDave May PetscFunctionReturn(0); 252d8588912SDave May } 253d8588912SDave May 254d8588912SDave May #undef __FUNCT__ 255f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 256f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 257d8588912SDave May { 258207556f9SJed Brown PetscErrorCode ierr; 259f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 260f349c1fdSJed Brown PetscInt j; 261f349c1fdSJed Brown Mat sub; 262d8588912SDave May 263d8588912SDave May PetscFunctionBegin; 2641db2f0e5SJed Brown sub = (row < vs->nc) ? vs->m[row][row] : (Mat)PETSC_NULL; /* Prefer to find on the diagonal */ 265f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 2664994cf47SJed Brown if (sub) {ierr = MatSetUp(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__ "MatNestFindNonzeroSubMatCol" 273f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 274f349c1fdSJed Brown { 275207556f9SJed Brown PetscErrorCode ierr; 276f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 277f349c1fdSJed Brown PetscInt i; 278f349c1fdSJed Brown Mat sub; 279f349c1fdSJed Brown 280f349c1fdSJed Brown PetscFunctionBegin; 2811db2f0e5SJed Brown sub = (col < vs->nr) ? vs->m[col][col] : (Mat)PETSC_NULL; /* Prefer to find on the diagonal */ 282f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 2834994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 284f349c1fdSJed Brown *B = sub; 285f349c1fdSJed Brown PetscFunctionReturn(0); 286d8588912SDave May } 287d8588912SDave May 288f349c1fdSJed Brown #undef __FUNCT__ 289f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 290f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 291f349c1fdSJed Brown { 292f349c1fdSJed Brown PetscErrorCode ierr; 293f349c1fdSJed Brown PetscInt i; 294f349c1fdSJed Brown PetscBool flg; 295f349c1fdSJed Brown 296f349c1fdSJed Brown PetscFunctionBegin; 297f349c1fdSJed Brown PetscValidPointer(list,3); 298f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 299f349c1fdSJed Brown PetscValidIntPointer(found,5); 300f349c1fdSJed Brown *found = -1; 301f349c1fdSJed Brown for (i=0; i<n; i++) { 302207556f9SJed Brown if (!list[i]) continue; 303f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 304f349c1fdSJed Brown if (flg) { 305f349c1fdSJed Brown *found = i; 306f349c1fdSJed Brown PetscFunctionReturn(0); 307f349c1fdSJed Brown } 308f349c1fdSJed Brown } 309f349c1fdSJed Brown SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set"); 310f349c1fdSJed Brown PetscFunctionReturn(0); 311f349c1fdSJed Brown } 312f349c1fdSJed Brown 313f349c1fdSJed Brown #undef __FUNCT__ 3148188e55aSJed Brown #define __FUNCT__ "MatNestGetRow" 3158188e55aSJed Brown /* Get a block row as a new MatNest */ 3168188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 3178188e55aSJed Brown { 3188188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3198188e55aSJed Brown char keyname[256]; 3208188e55aSJed Brown PetscErrorCode ierr; 3218188e55aSJed Brown 3228188e55aSJed Brown PetscFunctionBegin; 3238188e55aSJed Brown *B = PETSC_NULL; 3248caf3d72SBarry Smith ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr); 3258188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 3268188e55aSJed Brown if (*B) PetscFunctionReturn(0); 3278188e55aSJed Brown 3288188e55aSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 3298188e55aSJed Brown (*B)->assembled = A->assembled; 3308188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 3318188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 3328188e55aSJed Brown PetscFunctionReturn(0); 3338188e55aSJed Brown } 3348188e55aSJed Brown 3358188e55aSJed Brown #undef __FUNCT__ 336f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 337f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 338f349c1fdSJed Brown { 339f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3408188e55aSJed Brown PetscErrorCode ierr; 3416b3a5b13SJed Brown PetscInt row,col; 342e072481dSJed Brown PetscBool same,isFullCol,isFullColGlobal; 343f349c1fdSJed Brown 344f349c1fdSJed Brown PetscFunctionBegin; 3458188e55aSJed Brown /* Check if full column space. This is a hack */ 3468188e55aSJed Brown isFullCol = PETSC_FALSE; 347251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 3488188e55aSJed Brown if (same) { 34977019fcaSJed Brown PetscInt n,first,step,i,an,am,afirst,astep; 3508188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 3518188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 35277019fcaSJed Brown isFullCol = PETSC_TRUE; 35305ce4453SJed Brown for (i=0,an=A->cmap->rstart; i<vs->nc; i++) { 35477019fcaSJed Brown ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr); 35577019fcaSJed Brown ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr); 35677019fcaSJed Brown if (afirst != an || astep != step) isFullCol = PETSC_FALSE; 35777019fcaSJed Brown an += am; 35877019fcaSJed Brown } 35905ce4453SJed Brown if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE; 3608188e55aSJed Brown } 361c3aae356SJed Brown ierr = MPI_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,((PetscObject)iscol)->comm);CHKERRQ(ierr); 3628188e55aSJed Brown 363e072481dSJed Brown if (isFullColGlobal) { 3648188e55aSJed Brown PetscInt row; 3658188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 3668188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 3678188e55aSJed Brown } else { 368f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 369f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 370f349c1fdSJed Brown *B = vs->m[row][col]; 3718188e55aSJed Brown } 372f349c1fdSJed Brown PetscFunctionReturn(0); 373f349c1fdSJed Brown } 374f349c1fdSJed Brown 375f349c1fdSJed Brown #undef __FUNCT__ 376f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 377f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 378f349c1fdSJed Brown { 379f349c1fdSJed Brown PetscErrorCode ierr; 380f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 381f349c1fdSJed Brown Mat sub; 382f349c1fdSJed Brown 383f349c1fdSJed Brown PetscFunctionBegin; 384f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 385f349c1fdSJed Brown switch (reuse) { 386f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 3877874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 388f349c1fdSJed Brown *B = sub; 389f349c1fdSJed Brown break; 390f349c1fdSJed Brown case MAT_REUSE_MATRIX: 391f349c1fdSJed Brown if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 392f349c1fdSJed Brown break; 393f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 394f349c1fdSJed Brown break; 395f349c1fdSJed Brown } 396f349c1fdSJed Brown PetscFunctionReturn(0); 397f349c1fdSJed Brown } 398f349c1fdSJed Brown 399f349c1fdSJed Brown #undef __FUNCT__ 400f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 401f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 402f349c1fdSJed Brown { 403f349c1fdSJed Brown PetscErrorCode ierr; 404f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 405f349c1fdSJed Brown Mat sub; 406f349c1fdSJed Brown 407f349c1fdSJed Brown PetscFunctionBegin; 408f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 409f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 410f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 411f349c1fdSJed Brown *B = sub; 412d8588912SDave May PetscFunctionReturn(0); 413d8588912SDave May } 414d8588912SDave May 415d8588912SDave May #undef __FUNCT__ 416d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 417207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 418d8588912SDave May { 419d8588912SDave May PetscErrorCode ierr; 420f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 421f349c1fdSJed Brown Mat sub; 422d8588912SDave May 423d8588912SDave May PetscFunctionBegin; 424f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 425f349c1fdSJed Brown if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 426f349c1fdSJed Brown if (sub) { 427f349c1fdSJed Brown if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4286bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 429d8588912SDave May } 430d8588912SDave May PetscFunctionReturn(0); 431d8588912SDave May } 432d8588912SDave May 433d8588912SDave May #undef __FUNCT__ 4347874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 4357874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 4367874fa86SDave May { 4377874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 4387874fa86SDave May PetscInt i; 4397874fa86SDave May PetscErrorCode ierr; 4407874fa86SDave May 4417874fa86SDave May PetscFunctionBegin; 4427874fa86SDave May for (i=0; i<bA->nr; i++) { 443429bac76SJed Brown Vec bv; 444429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4457874fa86SDave May if (bA->m[i][i]) { 446429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 4477874fa86SDave May } else { 448429bac76SJed Brown ierr = VecSet(bv,1.0);CHKERRQ(ierr); 4497874fa86SDave May } 450429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 4517874fa86SDave May } 4527874fa86SDave May PetscFunctionReturn(0); 4537874fa86SDave May } 4547874fa86SDave May 4557874fa86SDave May #undef __FUNCT__ 4567874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 4577874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 4587874fa86SDave May { 4597874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 460429bac76SJed Brown Vec bl,*br; 4617874fa86SDave May PetscInt i,j; 4627874fa86SDave May PetscErrorCode ierr; 4637874fa86SDave May 4647874fa86SDave May PetscFunctionBegin; 465429bac76SJed Brown ierr = PetscMalloc(bA->nc*sizeof(Vec),&br);CHKERRQ(ierr); 466429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 4677874fa86SDave May for (i=0; i<bA->nr; i++) { 468429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4697874fa86SDave May for (j=0; j<bA->nc; j++) { 4707874fa86SDave May if (bA->m[i][j]) { 471429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 4727874fa86SDave May } 4737874fa86SDave May } 474a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 4757874fa86SDave May } 476429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 477429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 4787874fa86SDave May PetscFunctionReturn(0); 4797874fa86SDave May } 4807874fa86SDave May 4817874fa86SDave May #undef __FUNCT__ 482a061e289SJed Brown #define __FUNCT__ "MatScale_Nest" 483a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 484a061e289SJed Brown { 485a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 486a061e289SJed Brown PetscInt i,j; 487a061e289SJed Brown PetscErrorCode ierr; 488a061e289SJed Brown 489a061e289SJed Brown PetscFunctionBegin; 490a061e289SJed Brown for (i=0; i<bA->nr; i++) { 491a061e289SJed Brown for (j=0; j<bA->nc; j++) { 492a061e289SJed Brown if (bA->m[i][j]) { 493a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 494a061e289SJed Brown } 495a061e289SJed Brown } 496a061e289SJed Brown } 497a061e289SJed Brown PetscFunctionReturn(0); 498a061e289SJed Brown } 499a061e289SJed Brown 500a061e289SJed Brown #undef __FUNCT__ 501a061e289SJed Brown #define __FUNCT__ "MatShift_Nest" 502a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 503a061e289SJed Brown { 504a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 505a061e289SJed Brown PetscInt i; 506a061e289SJed Brown PetscErrorCode ierr; 507a061e289SJed Brown 508a061e289SJed Brown PetscFunctionBegin; 509a061e289SJed Brown for (i=0; i<bA->nr; i++) { 510a061e289SJed 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); 511a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 512a061e289SJed Brown } 513a061e289SJed Brown PetscFunctionReturn(0); 514a061e289SJed Brown } 515a061e289SJed Brown 516a061e289SJed Brown #undef __FUNCT__ 517d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest" 518207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left) 519d8588912SDave May { 520d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 521d8588912SDave May Vec *L,*R; 522d8588912SDave May MPI_Comm comm; 523d8588912SDave May PetscInt i,j; 524d8588912SDave May PetscErrorCode ierr; 525d8588912SDave May 526d8588912SDave May PetscFunctionBegin; 527d8588912SDave May comm = ((PetscObject)A)->comm; 528d8588912SDave May if (right) { 529d8588912SDave May /* allocate R */ 530d8588912SDave May ierr = PetscMalloc(sizeof(Vec) * bA->nc, &R);CHKERRQ(ierr); 531d8588912SDave May /* Create the right vectors */ 532d8588912SDave May for (j=0; j<bA->nc; j++) { 533d8588912SDave May for (i=0; i<bA->nr; i++) { 534d8588912SDave May if (bA->m[i][j]) { 535d8588912SDave May ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr); 536d8588912SDave May break; 537d8588912SDave May } 538d8588912SDave May } 539d8588912SDave May if (i==bA->nr) { 540d8588912SDave May /* have an empty column */ 541d8588912SDave May SETERRQ(((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 542d8588912SDave May } 543d8588912SDave May } 544f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 545d8588912SDave May /* hand back control to the nest vector */ 546d8588912SDave May for (j=0; j<bA->nc; j++) { 5476bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 548d8588912SDave May } 549d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 550d8588912SDave May } 551d8588912SDave May 552d8588912SDave May if (left) { 553d8588912SDave May /* allocate L */ 554d8588912SDave May ierr = PetscMalloc(sizeof(Vec) * bA->nr, &L);CHKERRQ(ierr); 555d8588912SDave May /* Create the left vectors */ 556d8588912SDave May for (i=0; i<bA->nr; i++) { 557d8588912SDave May for (j=0; j<bA->nc; j++) { 558d8588912SDave May if (bA->m[i][j]) { 559d8588912SDave May ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr); 560d8588912SDave May break; 561d8588912SDave May } 562d8588912SDave May } 563d8588912SDave May if (j==bA->nc) { 564d8588912SDave May /* have an empty row */ 565d8588912SDave May SETERRQ(((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 566d8588912SDave May } 567d8588912SDave May } 568d8588912SDave May 569f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 570d8588912SDave May for (i=0; i<bA->nr; i++) { 5716bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 572d8588912SDave May } 573d8588912SDave May 574d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 575d8588912SDave May } 576d8588912SDave May PetscFunctionReturn(0); 577d8588912SDave May } 578d8588912SDave May 579d8588912SDave May #undef __FUNCT__ 580d8588912SDave May #define __FUNCT__ "MatView_Nest" 581207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 582d8588912SDave May { 583d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 584d8588912SDave May PetscBool isascii; 585d8588912SDave May PetscInt i,j; 586d8588912SDave May PetscErrorCode ierr; 587d8588912SDave May 588d8588912SDave May PetscFunctionBegin; 589251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 590d8588912SDave May if (isascii) { 591d8588912SDave May 592d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n"); 593d8588912SDave May PetscViewerASCIIPushTab(viewer); /* push0 */ 594d8588912SDave May PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 595d8588912SDave May 596d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n"); 597d8588912SDave May for (i=0; i<bA->nr; i++) { 598d8588912SDave May for (j=0; j<bA->nc; j++) { 59919fd82e9SBarry Smith MatType type; 600270f95d7SJed Brown char name[256] = "",prefix[256] = ""; 601d8588912SDave May PetscInt NR,NC; 602d8588912SDave May PetscBool isNest = PETSC_FALSE; 603d8588912SDave May 604d8588912SDave May if (!bA->m[i][j]) { 605d8588912SDave May PetscViewerASCIIPrintf(viewer, "(%D,%D) : PETSC_NULL \n",i,j); 606d8588912SDave May continue; 607d8588912SDave May } 608d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 609d8588912SDave May ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr); 6108caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 6118caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 612251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 613d8588912SDave May 614270f95d7SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 615d8588912SDave May 616d8588912SDave May if (isNest) { 617270f95d7SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 618d8588912SDave May ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 619270f95d7SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 620d8588912SDave May } 621d8588912SDave May } 622d8588912SDave May } 623d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 624d8588912SDave May } 625d8588912SDave May PetscFunctionReturn(0); 626d8588912SDave May } 627d8588912SDave May 628d8588912SDave May #undef __FUNCT__ 629d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 630207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 631d8588912SDave May { 632d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 633d8588912SDave May PetscInt i,j; 634d8588912SDave May PetscErrorCode ierr; 635d8588912SDave May 636d8588912SDave May PetscFunctionBegin; 637d8588912SDave May for (i=0; i<bA->nr; i++) { 638d8588912SDave May for (j=0; j<bA->nc; j++) { 639d8588912SDave May if (!bA->m[i][j]) continue; 640d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 641d8588912SDave May } 642d8588912SDave May } 643d8588912SDave May PetscFunctionReturn(0); 644d8588912SDave May } 645d8588912SDave May 646d8588912SDave May #undef __FUNCT__ 647d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 648207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 649d8588912SDave May { 650d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 651841e96a3SJed Brown Mat *b; 652841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 653d8588912SDave May PetscErrorCode ierr; 654d8588912SDave May 655d8588912SDave May PetscFunctionBegin; 656841e96a3SJed Brown ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr); 657841e96a3SJed Brown for (i=0; i<nr; i++) { 658841e96a3SJed Brown for (j=0; j<nc; j++) { 659841e96a3SJed Brown if (bA->m[i][j]) { 660841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 661841e96a3SJed Brown } else { 662841e96a3SJed Brown b[i*nc+j] = PETSC_NULL; 663d8588912SDave May } 664d8588912SDave May } 665d8588912SDave May } 666f349c1fdSJed Brown ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 667841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 668841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 6696bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 670d8588912SDave May } 671d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 672d8588912SDave May 673841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 674841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 675d8588912SDave May PetscFunctionReturn(0); 676d8588912SDave May } 677d8588912SDave May 678d8588912SDave May /* nest api */ 679d8588912SDave May EXTERN_C_BEGIN 680d8588912SDave May #undef __FUNCT__ 681d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 682d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 683d8588912SDave May { 684d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 6855fd66863SKarl Rupp 686d8588912SDave May PetscFunctionBegin; 687d8588912SDave May if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 688d8588912SDave May if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 689d8588912SDave May *mat = bA->m[idxm][jdxm]; 690d8588912SDave May PetscFunctionReturn(0); 691d8588912SDave May } 692d8588912SDave May EXTERN_C_END 693d8588912SDave May 694d8588912SDave May #undef __FUNCT__ 695d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 6969ba0d327SJed Brown /*@ 697d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 698d8588912SDave May 699d8588912SDave May Not collective 700d8588912SDave May 701d8588912SDave May Input Parameters: 702629881c0SJed Brown + A - nest matrix 703d8588912SDave May . idxm - index of the matrix within the nest matrix 704629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 705d8588912SDave May 706d8588912SDave May Output Parameter: 707d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 708d8588912SDave May 709d8588912SDave May Level: developer 710d8588912SDave May 711d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 712d8588912SDave May @*/ 7137087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 714d8588912SDave May { 715699a902aSJed Brown PetscErrorCode ierr; 716d8588912SDave May 717d8588912SDave May PetscFunctionBegin; 718699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 719d8588912SDave May PetscFunctionReturn(0); 720d8588912SDave May } 721d8588912SDave May 722d8588912SDave May EXTERN_C_BEGIN 723d8588912SDave May #undef __FUNCT__ 7240782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 7250782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 7260782ca92SJed Brown { 7270782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 7280782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 7290782ca92SJed Brown PetscErrorCode ierr; 7300782ca92SJed Brown 7310782ca92SJed Brown PetscFunctionBegin; 7320782ca92SJed Brown if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 7330782ca92SJed Brown if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 7340782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 7350782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 7360782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 7370782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 7380782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 7390782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 7400782ca92SJed 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); 7410782ca92SJed 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); 7420782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 7430782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 7440782ca92SJed Brown bA->m[idxm][jdxm] = mat; 7450782ca92SJed Brown PetscFunctionReturn(0); 7460782ca92SJed Brown } 7470782ca92SJed Brown EXTERN_C_END 7480782ca92SJed Brown 7490782ca92SJed Brown #undef __FUNCT__ 7500782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 7519ba0d327SJed Brown /*@ 7520782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 7530782ca92SJed Brown 7540782ca92SJed Brown Logically collective on the submatrix communicator 7550782ca92SJed Brown 7560782ca92SJed Brown Input Parameters: 7570782ca92SJed Brown + A - nest matrix 7580782ca92SJed Brown . idxm - index of the matrix within the nest matrix 7590782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 7600782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 7610782ca92SJed Brown 7620782ca92SJed Brown Notes: 7630782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 7640782ca92SJed Brown 7650782ca92SJed Brown This increments the reference count of the submatrix. 7660782ca92SJed Brown 7670782ca92SJed Brown Level: developer 7680782ca92SJed Brown 7690782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 7700782ca92SJed Brown @*/ 7710782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 7720782ca92SJed Brown { 7730782ca92SJed Brown PetscErrorCode ierr; 7740782ca92SJed Brown 7750782ca92SJed Brown PetscFunctionBegin; 7760782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 7770782ca92SJed Brown PetscFunctionReturn(0); 7780782ca92SJed Brown } 7790782ca92SJed Brown 7800782ca92SJed Brown EXTERN_C_BEGIN 7810782ca92SJed Brown #undef __FUNCT__ 782d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 783d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 784d8588912SDave May { 785d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 7865fd66863SKarl Rupp 787d8588912SDave May PetscFunctionBegin; 788d8588912SDave May if (M) { *M = bA->nr; } 789d8588912SDave May if (N) { *N = bA->nc; } 790d8588912SDave May if (mat) { *mat = bA->m; } 791d8588912SDave May PetscFunctionReturn(0); 792d8588912SDave May } 793d8588912SDave May EXTERN_C_END 794d8588912SDave May 795d8588912SDave May #undef __FUNCT__ 796d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 797d8588912SDave May /*@C 798d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 799d8588912SDave May 800d8588912SDave May Not collective 801d8588912SDave May 802d8588912SDave May Input Parameters: 803629881c0SJed Brown . A - nest matrix 804d8588912SDave May 805d8588912SDave May Output Parameter: 806629881c0SJed Brown + M - number of rows in the nest matrix 807d8588912SDave May . N - number of cols in the nest matrix 808629881c0SJed Brown - mat - 2d array of matrices 809d8588912SDave May 810d8588912SDave May Notes: 811d8588912SDave May 812d8588912SDave May The user should not free the array mat. 813d8588912SDave May 814d8588912SDave May Level: developer 815d8588912SDave May 816d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 817d8588912SDave May @*/ 8187087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 819d8588912SDave May { 820699a902aSJed Brown PetscErrorCode ierr; 821d8588912SDave May 822d8588912SDave May PetscFunctionBegin; 823699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 824d8588912SDave May PetscFunctionReturn(0); 825d8588912SDave May } 826d8588912SDave May 827d8588912SDave May EXTERN_C_BEGIN 828d8588912SDave May #undef __FUNCT__ 829d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 8307087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 831d8588912SDave May { 832d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 833d8588912SDave May 834d8588912SDave May PetscFunctionBegin; 835d8588912SDave May if (M) { *M = bA->nr; } 836d8588912SDave May if (N) { *N = bA->nc; } 837d8588912SDave May PetscFunctionReturn(0); 838d8588912SDave May } 839d8588912SDave May EXTERN_C_END 840d8588912SDave May 841d8588912SDave May #undef __FUNCT__ 842d8588912SDave May #define __FUNCT__ "MatNestGetSize" 8439ba0d327SJed Brown /*@ 844d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 845d8588912SDave May 846d8588912SDave May Not collective 847d8588912SDave May 848d8588912SDave May Input Parameters: 849d8588912SDave May . A - nest matrix 850d8588912SDave May 851d8588912SDave May Output Parameter: 852629881c0SJed Brown + M - number of rows in the nested mat 853629881c0SJed Brown - N - number of cols in the nested mat 854d8588912SDave May 855d8588912SDave May Notes: 856d8588912SDave May 857d8588912SDave May Level: developer 858d8588912SDave May 859d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 860d8588912SDave May @*/ 8617087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 862d8588912SDave May { 863699a902aSJed Brown PetscErrorCode ierr; 864d8588912SDave May 865d8588912SDave May PetscFunctionBegin; 866699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 867d8588912SDave May PetscFunctionReturn(0); 868d8588912SDave May } 869d8588912SDave May 870900e7ff2SJed Brown #undef __FUNCT__ 871900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest" 872900e7ff2SJed Brown PETSC_EXTERN_C PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 873900e7ff2SJed Brown { 874900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 875900e7ff2SJed Brown PetscInt i; 876900e7ff2SJed Brown 877900e7ff2SJed Brown PetscFunctionBegin; 878900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 879900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 880900e7ff2SJed Brown PetscFunctionReturn(0); 881900e7ff2SJed Brown } 882900e7ff2SJed Brown 883900e7ff2SJed Brown #undef __FUNCT__ 884900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs" 885900e7ff2SJed Brown /*@C 886900e7ff2SJed Brown MatNestGetISs - Returns the index sets partitioning the row and column spaces 887900e7ff2SJed Brown 888900e7ff2SJed Brown Not collective 889900e7ff2SJed Brown 890900e7ff2SJed Brown Input Parameters: 891900e7ff2SJed Brown . A - nest matrix 892900e7ff2SJed Brown 893900e7ff2SJed Brown Output Parameter: 894900e7ff2SJed Brown + rows - array of row index sets 895900e7ff2SJed Brown - cols - array of column index sets 896900e7ff2SJed Brown 897900e7ff2SJed Brown Level: advanced 898900e7ff2SJed Brown 899900e7ff2SJed Brown Notes: 900900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 901900e7ff2SJed Brown 902900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs() 903900e7ff2SJed Brown @*/ 904900e7ff2SJed Brown PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 905900e7ff2SJed Brown { 906900e7ff2SJed Brown PetscErrorCode ierr; 907900e7ff2SJed Brown 908900e7ff2SJed Brown PetscFunctionBegin; 909900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 910900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 911900e7ff2SJed Brown PetscFunctionReturn(0); 912900e7ff2SJed Brown } 913900e7ff2SJed Brown 914900e7ff2SJed Brown #undef __FUNCT__ 915900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest" 916900e7ff2SJed Brown PETSC_EXTERN_C PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 917900e7ff2SJed Brown { 918900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 919900e7ff2SJed Brown PetscInt i; 920900e7ff2SJed Brown 921900e7ff2SJed Brown PetscFunctionBegin; 922900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 923900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 924900e7ff2SJed Brown PetscFunctionReturn(0); 925900e7ff2SJed Brown } 926900e7ff2SJed Brown 927900e7ff2SJed Brown #undef __FUNCT__ 928900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs" 929900e7ff2SJed Brown /*@C 930900e7ff2SJed Brown MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 931900e7ff2SJed Brown 932900e7ff2SJed Brown Not collective 933900e7ff2SJed Brown 934900e7ff2SJed Brown Input Parameters: 935900e7ff2SJed Brown . A - nest matrix 936900e7ff2SJed Brown 937900e7ff2SJed Brown Output Parameter: 938900e7ff2SJed Brown + rows - array of row index sets (or PETSC_NULL to ignore) 939900e7ff2SJed Brown - cols - array of column index sets (or PETSC_NULL to ignore) 940900e7ff2SJed Brown 941900e7ff2SJed Brown Level: advanced 942900e7ff2SJed Brown 943900e7ff2SJed Brown Notes: 944900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 945900e7ff2SJed Brown 946900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs() 947900e7ff2SJed Brown @*/ 948900e7ff2SJed Brown PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 949900e7ff2SJed Brown { 950900e7ff2SJed Brown PetscErrorCode ierr; 951900e7ff2SJed Brown 952900e7ff2SJed Brown PetscFunctionBegin; 953900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 954900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 955900e7ff2SJed Brown PetscFunctionReturn(0); 956900e7ff2SJed Brown } 957900e7ff2SJed Brown 958207556f9SJed Brown EXTERN_C_BEGIN 959207556f9SJed Brown #undef __FUNCT__ 960207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 96119fd82e9SBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,VecType vtype) 962207556f9SJed Brown { 963207556f9SJed Brown PetscErrorCode ierr; 964207556f9SJed Brown PetscBool flg; 965207556f9SJed Brown 966207556f9SJed Brown PetscFunctionBegin; 967207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 968207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 96912b53f24SSatish Balay if (flg) A->ops->getvecs = MatGetVecs_Nest; 97012b53f24SSatish Balay else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*))0; 971207556f9SJed Brown PetscFunctionReturn(0); 972207556f9SJed Brown } 973207556f9SJed Brown EXTERN_C_END 974207556f9SJed Brown 975207556f9SJed Brown #undef __FUNCT__ 976207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 977207556f9SJed Brown /*@C 978207556f9SJed Brown MatNestSetVecType - Sets the type of Vec returned by MatGetVecs() 979207556f9SJed Brown 980207556f9SJed Brown Not collective 981207556f9SJed Brown 982207556f9SJed Brown Input Parameters: 983207556f9SJed Brown + A - nest matrix 984207556f9SJed Brown - vtype - type to use for creating vectors 985207556f9SJed Brown 986207556f9SJed Brown Notes: 987207556f9SJed Brown 988207556f9SJed Brown Level: developer 989207556f9SJed Brown 990207556f9SJed Brown .seealso: MatGetVecs() 991207556f9SJed Brown @*/ 99219fd82e9SBarry Smith PetscErrorCode MatNestSetVecType(Mat A,VecType vtype) 993207556f9SJed Brown { 994207556f9SJed Brown PetscErrorCode ierr; 995207556f9SJed Brown 996207556f9SJed Brown PetscFunctionBegin; 99719fd82e9SBarry Smith ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr); 998207556f9SJed Brown PetscFunctionReturn(0); 999207556f9SJed Brown } 1000207556f9SJed Brown 1001c8883902SJed Brown EXTERN_C_BEGIN 1002d8588912SDave May #undef __FUNCT__ 1003c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 1004c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1005d8588912SDave May { 1006c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 1007c8883902SJed Brown PetscInt i,j,m,n,M,N; 1008d8588912SDave May PetscErrorCode ierr; 1009d8588912SDave May 1010d8588912SDave May PetscFunctionBegin; 1011c8883902SJed Brown s->nr = nr; 1012c8883902SJed Brown s->nc = nc; 1013d8588912SDave May 1014c8883902SJed Brown /* Create space for submatrices */ 1015c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr); 1016c8883902SJed Brown for (i=0; i<nr; i++) { 1017c8883902SJed Brown ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr); 1018d8588912SDave May } 1019c8883902SJed Brown for (i=0; i<nr; i++) { 1020c8883902SJed Brown for (j=0; j<nc; j++) { 1021c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1022c8883902SJed Brown if (a[i*nc+j]) { 1023c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1024d8588912SDave May } 1025d8588912SDave May } 1026d8588912SDave May } 1027d8588912SDave May 10288188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1029d8588912SDave May 1030c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr); 1031c8883902SJed Brown ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr); 1032c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1033c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1034d8588912SDave May 10358188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1036d8588912SDave May 1037c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1038c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1039c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1040c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1041c8883902SJed Brown 1042c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1043c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1044c8883902SJed Brown 1045c8883902SJed Brown ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr); 1046c8883902SJed Brown ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr); 1047c8883902SJed Brown ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr); 1048d8588912SDave May PetscFunctionReturn(0); 1049d8588912SDave May } 1050c8883902SJed Brown EXTERN_C_END 1051d8588912SDave May 1052c8883902SJed Brown #undef __FUNCT__ 1053c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 1054c8883902SJed Brown /*@ 1055c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1056c8883902SJed Brown 1057c8883902SJed Brown Collective on Mat 1058c8883902SJed Brown 1059c8883902SJed Brown Input Parameter: 1060c8883902SJed Brown + N - nested matrix 1061c8883902SJed Brown . nr - number of nested row blocks 1062c8883902SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1063c8883902SJed Brown . nc - number of nested column blocks 1064c8883902SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1065c8883902SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1066c8883902SJed Brown 1067c8883902SJed Brown Level: advanced 1068c8883902SJed Brown 1069c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 1070c8883902SJed Brown @*/ 1071c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1072c8883902SJed Brown { 1073c8883902SJed Brown PetscErrorCode ierr; 1074c8883902SJed Brown PetscInt i; 1075c8883902SJed Brown 1076c8883902SJed Brown PetscFunctionBegin; 1077c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1078c8883902SJed Brown if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1079c8883902SJed Brown if (nr && is_row) { 1080c8883902SJed Brown PetscValidPointer(is_row,3); 1081c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1082c8883902SJed Brown } 1083c8883902SJed Brown if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 10841664e352SJed Brown if (nc && is_col) { 1085c8883902SJed Brown PetscValidPointer(is_col,5); 1086c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1087c8883902SJed Brown } 1088c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 1089c8883902SJed 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); 1090c8883902SJed Brown PetscFunctionReturn(0); 1091c8883902SJed Brown } 1092d8588912SDave May 109377019fcaSJed Brown #undef __FUNCT__ 109477019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private" 109577019fcaSJed Brown static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb) 109677019fcaSJed Brown { 109777019fcaSJed Brown PetscErrorCode ierr; 109877019fcaSJed Brown PetscBool flg; 109977019fcaSJed Brown PetscInt i,j,m,mi,*ix; 110077019fcaSJed Brown 110177019fcaSJed Brown PetscFunctionBegin; 110277019fcaSJed Brown for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 110377019fcaSJed Brown if (islocal[i]) { 110477019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 110577019fcaSJed Brown flg = PETSC_TRUE; /* We found a non-trivial entry */ 110677019fcaSJed Brown } else { 110777019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 110877019fcaSJed Brown } 110977019fcaSJed Brown m += mi; 111077019fcaSJed Brown } 111177019fcaSJed Brown if (flg) { 111277019fcaSJed Brown ierr = PetscMalloc(m*sizeof(*ix),&ix);CHKERRQ(ierr); 111377019fcaSJed Brown for (i=0,n=0; i<n; i++) { 111477019fcaSJed Brown ISLocalToGlobalMapping smap = PETSC_NULL; 111577019fcaSJed Brown VecScatter scat; 111677019fcaSJed Brown IS isreq; 111777019fcaSJed Brown Vec lvec,gvec; 11183361c9a7SJed Brown union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x; 111977019fcaSJed Brown Mat sub; 112077019fcaSJed Brown 112177019fcaSJed Brown if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support when scalars smaller than integers"); 112277019fcaSJed Brown if (colflg) { 112377019fcaSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 112477019fcaSJed Brown } else { 112577019fcaSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 112677019fcaSJed Brown } 112777019fcaSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,PETSC_NULL);CHKERRQ(ierr);} 112877019fcaSJed Brown if (islocal[i]) { 112977019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 113077019fcaSJed Brown } else { 113177019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 113277019fcaSJed Brown } 113377019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = j; 113477019fcaSJed Brown if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);} 113577019fcaSJed Brown /* 113677019fcaSJed Brown Now we need to extract the monolithic global indices that correspond to the given split global indices. 113777019fcaSJed 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. 113877019fcaSJed Brown The approach here is ugly because it uses VecScatter to move indices. 113977019fcaSJed Brown */ 114077019fcaSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr); 114177019fcaSJed Brown ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr); 114277019fcaSJed Brown ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr); 114377019fcaSJed Brown ierr = VecScatterCreate(gvec,isreq,lvec,PETSC_NULL,&scat);CHKERRQ(ierr); 114477019fcaSJed Brown ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 114577019fcaSJed Brown for (j=0; j<mi; j++) x[j].integer = ix[m+j]; 114677019fcaSJed Brown ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 114777019fcaSJed Brown ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 114877019fcaSJed Brown ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 114977019fcaSJed Brown ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 115077019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = x[j].integer; 115177019fcaSJed Brown ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 115277019fcaSJed Brown ierr = VecDestroy(&lvec);CHKERRQ(ierr); 115377019fcaSJed Brown ierr = VecDestroy(&gvec);CHKERRQ(ierr); 115477019fcaSJed Brown ierr = ISDestroy(&isreq);CHKERRQ(ierr); 115577019fcaSJed Brown ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 115677019fcaSJed Brown m += mi; 115777019fcaSJed Brown } 115877019fcaSJed Brown ierr = ISLocalToGlobalMappingCreate(((PetscObject)A)->comm,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 115977019fcaSJed Brown *ltogb = PETSC_NULL; 116077019fcaSJed Brown } else { 116177019fcaSJed Brown *ltog = PETSC_NULL; 116277019fcaSJed Brown *ltogb = PETSC_NULL; 116377019fcaSJed Brown } 116477019fcaSJed Brown PetscFunctionReturn(0); 116577019fcaSJed Brown } 116677019fcaSJed Brown 116777019fcaSJed Brown 1168d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1169d8588912SDave May /* 1170d8588912SDave May nprocessors = NP 1171d8588912SDave May Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1)) 1172d8588912SDave May proc 0: => (g_0,h_0,) 1173d8588912SDave May proc 1: => (g_1,h_1,) 1174d8588912SDave May ... 1175d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1176d8588912SDave May 1177d8588912SDave May proc 0: proc 1: proc nprocs-1: 1178d8588912SDave May is[0] = (0,1,2,...,nlocal(g_0)-1) (0,1,...,nlocal(g_1)-1) (0,1,...,nlocal(g_NP-1)) 1179d8588912SDave May 1180d8588912SDave May proc 0: 1181d8588912SDave May is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1) 1182d8588912SDave May proc 1: 1183d8588912SDave May is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1) 1184d8588912SDave May 1185d8588912SDave May proc NP-1: 1186d8588912SDave May is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1) 1187d8588912SDave May */ 1188d8588912SDave May #undef __FUNCT__ 1189d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1190841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1191d8588912SDave May { 1192e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 11938188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1194d8588912SDave May PetscErrorCode ierr; 119587743a01SMatthew G Knepley Mat sub = PETSC_NULL; 1196d8588912SDave May 1197d8588912SDave May PetscFunctionBegin; 11988188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr); 11998188e55aSJed Brown ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr); 1200d8588912SDave May if (is_row) { /* valid IS is passed in */ 1201d8588912SDave May /* refs on is[] are incremeneted */ 1202e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1203d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 1204e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1205d8588912SDave May } 12062ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 12078188e55aSJed Brown nsum = 0; 12088188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 12098188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12108188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 12118188e55aSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 12128188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12138188e55aSJed Brown nsum += n; 12148188e55aSJed Brown } 121530bc264bSJed Brown ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 121630bc264bSJed Brown offset -= nsum; 1217e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1218f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12192ae74bdbSJed Brown ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr); 12202ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1221e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1222e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 12232ae74bdbSJed Brown offset += n; 1224d8588912SDave May } 1225d8588912SDave May } 1226d8588912SDave May 1227d8588912SDave May if (is_col) { /* valid IS is passed in */ 1228d8588912SDave May /* refs on is[] are incremeneted */ 1229e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1230d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 1231e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1232d8588912SDave May } 12332ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 12342ae74bdbSJed Brown offset = A->cmap->rstart; 12358188e55aSJed Brown nsum = 0; 12368188e55aSJed Brown for (j=0; j<vs->nc; j++) { 12378188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 12388188e55aSJed Brown if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 12398188e55aSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 12408188e55aSJed Brown if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 12418188e55aSJed Brown nsum += n; 12428188e55aSJed Brown } 124330bc264bSJed Brown ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 124430bc264bSJed Brown offset -= nsum; 1245e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1246f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 12472ae74bdbSJed Brown ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr); 12482ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1249e2d7f03fSJed Brown ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1250e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 12512ae74bdbSJed Brown offset += n; 1252d8588912SDave May } 1253d8588912SDave May } 1254e2d7f03fSJed Brown 1255e2d7f03fSJed Brown /* Set up the local ISs */ 1256e2d7f03fSJed Brown ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr); 1257e2d7f03fSJed Brown ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr); 1258e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1259e2d7f03fSJed Brown IS isloc; 12608188e55aSJed Brown ISLocalToGlobalMapping rmap = PETSC_NULL; 1261e2d7f03fSJed Brown PetscInt nlocal,bs; 1262e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 12638188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);} 1264207556f9SJed Brown if (rmap) { 1265e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1266e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1267e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1268e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1269207556f9SJed Brown } else { 1270207556f9SJed Brown nlocal = 0; 1271207556f9SJed Brown isloc = PETSC_NULL; 1272207556f9SJed Brown } 1273e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1274e2d7f03fSJed Brown offset += nlocal; 1275e2d7f03fSJed Brown } 12768188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1277e2d7f03fSJed Brown IS isloc; 12788188e55aSJed Brown ISLocalToGlobalMapping cmap = PETSC_NULL; 1279e2d7f03fSJed Brown PetscInt nlocal,bs; 1280e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 12818188e55aSJed Brown if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);} 1282207556f9SJed Brown if (cmap) { 1283e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1284e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1285e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1286e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1287207556f9SJed Brown } else { 1288207556f9SJed Brown nlocal = 0; 1289207556f9SJed Brown isloc = PETSC_NULL; 1290207556f9SJed Brown } 1291e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1292e2d7f03fSJed Brown offset += nlocal; 1293e2d7f03fSJed Brown } 12940189643fSJed Brown 129577019fcaSJed Brown /* Set up the aggregate ISLocalToGlobalMapping */ 129677019fcaSJed Brown { 129777019fcaSJed Brown ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb; 129877019fcaSJed Brown ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr); 129977019fcaSJed Brown ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr); 130077019fcaSJed Brown if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 130177019fcaSJed Brown if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);} 130277019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 130377019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr); 130477019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 130577019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr); 130677019fcaSJed Brown } 130777019fcaSJed Brown 13080189643fSJed Brown #if defined(PETSC_USE_DEBUG) 13090189643fSJed Brown for (i=0; i<vs->nr; i++) { 13100189643fSJed Brown for (j=0; j<vs->nc; j++) { 13110189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 13120189643fSJed Brown Mat B = vs->m[i][j]; 13130189643fSJed Brown if (!B) continue; 13140189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 13150189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 13160189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 13170189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 13180189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 13190189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 13200189643fSJed 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); 13210189643fSJed 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); 13220189643fSJed Brown } 13230189643fSJed Brown } 13240189643fSJed Brown #endif 1325a061e289SJed Brown 1326a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1327a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1328a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1329a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1330a061e289SJed Brown } 1331a061e289SJed Brown } 1332a061e289SJed Brown A->assembled = PETSC_TRUE; 1333d8588912SDave May PetscFunctionReturn(0); 1334d8588912SDave May } 1335d8588912SDave May 1336d8588912SDave May #undef __FUNCT__ 1337d8588912SDave May #define __FUNCT__ "MatCreateNest" 1338*45c38901SJed Brown /*@C 1339659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1340659c6bb0SJed Brown 1341659c6bb0SJed Brown Collective on Mat 1342659c6bb0SJed Brown 1343659c6bb0SJed Brown Input Parameter: 1344659c6bb0SJed Brown + comm - Communicator for the new Mat 1345659c6bb0SJed Brown . nr - number of nested row blocks 1346659c6bb0SJed Brown . is_row - index sets for each nested row block, or PETSC_NULL to make contiguous 1347659c6bb0SJed Brown . nc - number of nested column blocks 1348659c6bb0SJed Brown . is_col - index sets for each nested column block, or PETSC_NULL to make contiguous 1349659c6bb0SJed Brown - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL 1350659c6bb0SJed Brown 1351659c6bb0SJed Brown Output Parameter: 1352659c6bb0SJed Brown . B - new matrix 1353659c6bb0SJed Brown 1354659c6bb0SJed Brown Level: advanced 1355659c6bb0SJed Brown 1356950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST 1357659c6bb0SJed Brown @*/ 13587087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1359d8588912SDave May { 1360d8588912SDave May Mat A; 1361d8588912SDave May PetscErrorCode ierr; 1362d8588912SDave May 1363d8588912SDave May PetscFunctionBegin; 1364c8883902SJed Brown *B = 0; 1365d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1366c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 13677ae8954aSSatish Balay ierr = MatSetUp(A);CHKERRQ(ierr); 1368c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1369d8588912SDave May *B = A; 1370d8588912SDave May PetscFunctionReturn(0); 1371d8588912SDave May } 1372659c6bb0SJed Brown 1373629c3df2SDmitry Karpeev EXTERN_C_BEGIN 1374629c3df2SDmitry Karpeev #undef __FUNCT__ 1375629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ" 1376629c3df2SDmitry Karpeev PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1377629c3df2SDmitry Karpeev { 1378629c3df2SDmitry Karpeev PetscErrorCode ierr; 1379629c3df2SDmitry Karpeev Mat_Nest *nest = (Mat_Nest*)A->data; 1380629c3df2SDmitry Karpeev PetscInt m,n,M,N,i,j,k,*dnnz,*onnz; 1381629c3df2SDmitry Karpeev Mat C; 1382629c3df2SDmitry Karpeev 1383629c3df2SDmitry Karpeev PetscFunctionBegin; 1384629c3df2SDmitry Karpeev ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1385629c3df2SDmitry Karpeev ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1386629c3df2SDmitry Karpeev switch (reuse) { 1387629c3df2SDmitry Karpeev case MAT_INITIAL_MATRIX: 1388629c3df2SDmitry Karpeev ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 1389629c3df2SDmitry Karpeev ierr = MatSetType(C,newtype);CHKERRQ(ierr); 1390629c3df2SDmitry Karpeev ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 1391629c3df2SDmitry Karpeev *newmat = C; 1392629c3df2SDmitry Karpeev break; 1393629c3df2SDmitry Karpeev case MAT_REUSE_MATRIX: 1394629c3df2SDmitry Karpeev C = *newmat; 1395629c3df2SDmitry Karpeev break; 1396629c3df2SDmitry Karpeev default: SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatReuse"); 1397629c3df2SDmitry Karpeev } 1398629c3df2SDmitry Karpeev 1399629c3df2SDmitry Karpeev /* Preallocation */ 1400629c3df2SDmitry Karpeev ierr = PetscMalloc(2*m*sizeof(PetscInt),&dnnz);CHKERRQ(ierr); 1401629c3df2SDmitry Karpeev onnz = dnnz + m; 1402629c3df2SDmitry Karpeev for (k=0; k<m; k++) { 1403629c3df2SDmitry Karpeev dnnz[k] = 0; 1404629c3df2SDmitry Karpeev onnz[k] = 0; 1405629c3df2SDmitry Karpeev } 1406629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1407629c3df2SDmitry Karpeev IS bNis; 1408629c3df2SDmitry Karpeev PetscInt bN; 1409629c3df2SDmitry Karpeev const PetscInt *bNindices; 1410629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1411629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1412629c3df2SDmitry Karpeev ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr); 1413629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1414629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1415629c3df2SDmitry Karpeev PetscSF bmsf; 1416629c3df2SDmitry Karpeev PetscSFNode *bmedges; 1417629c3df2SDmitry Karpeev Mat B; 141834986ce2SMatthew G Knepley PetscInt bm, *bmdnnz, br; 1419629c3df2SDmitry Karpeev const PetscInt *bmindices; 1420629c3df2SDmitry Karpeev B = nest->m[i][j]; 1421629c3df2SDmitry Karpeev if (!B) continue; 1422629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1423629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1424629c3df2SDmitry Karpeev ierr = PetscSFCreate(((PetscObject)A)->comm, &bmsf);CHKERRQ(ierr); 1425629c3df2SDmitry Karpeev ierr = PetscMalloc(2*bm*sizeof(PetscSFNode),&bmedges);CHKERRQ(ierr); 1426629c3df2SDmitry Karpeev ierr = PetscMalloc(2*bm*sizeof(PetscInt),&bmdnnz);CHKERRQ(ierr); 1427629c3df2SDmitry Karpeev for (k = 0; k < 2*bm; ++k) bmdnnz[k] = 0; 1428629c3df2SDmitry Karpeev /* 1429629c3df2SDmitry Karpeev Locate the owners for all of the locally-owned global row indices for this row block. 1430629c3df2SDmitry Karpeev These determine the roots of PetscSF used to communicate preallocation data to row owners. 1431629c3df2SDmitry Karpeev The roots correspond to the dnnz and onnz entries; thus, there are two roots per row. 1432629c3df2SDmitry Karpeev */ 1433629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1434a4b3d3acSMatthew G Knepley PetscInt row = bmindices[br], rowowner = 0, brncols, col, colowner = 0; 1435629c3df2SDmitry Karpeev const PetscInt *brcols; 1436a4b3d3acSMatthew G Knepley PetscInt rowrel = 0;/* row's relative index on its owner rank */ 1437629c3df2SDmitry Karpeev PetscInt rowownerm; /* local row size on row's owning rank. */ 1438629c3df2SDmitry Karpeev ierr = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr); 1439629c3df2SDmitry Karpeev rowownerm = A->rmap->range[rowowner+1]-A->rmap->range[rowowner]; 1440629c3df2SDmitry Karpeev bmedges[br].rank = rowowner; bmedges[br].index = rowrel; /* edge from bmdnnz to dnnz */ 1441629c3df2SDmitry Karpeev bmedges[br].rank = rowowner; bmedges[br].index = rowrel+rowownerm; /* edge from bmonnz to onnz */ 1442629c3df2SDmitry Karpeev /* Now actually compute the data -- bmdnnz and bmonnz by looking at the global columns in the br row of this block. */ 1443629c3df2SDmitry Karpeev /* Note that this is not a pessimistic bound only because we assume the index sets embedding the blocks do not overlap. */ 1444629c3df2SDmitry Karpeev ierr = MatGetRow(B,br,&brncols,&brcols,PETSC_NULL);CHKERRQ(ierr); 1445629c3df2SDmitry Karpeev for (k=0; k<brncols; k++) { 1446629c3df2SDmitry Karpeev col = bNindices[brcols[k]]; 1447629c3df2SDmitry Karpeev ierr = PetscLayoutFindOwnerIndex(A->cmap,col,&colowner,PETSC_NULL);CHKERRQ(ierr); 1448629c3df2SDmitry Karpeev if (colowner == rowowner) bmdnnz[br]++; 1449629c3df2SDmitry Karpeev else onnz[br]++; 1450629c3df2SDmitry Karpeev } 1451629c3df2SDmitry Karpeev ierr = MatRestoreRow(B,br,&brncols,&brcols,PETSC_NULL);CHKERRQ(ierr); 1452629c3df2SDmitry Karpeev } 1453629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1454629c3df2SDmitry Karpeev /* bsf will have to take care of disposing of bedges. */ 1455629c3df2SDmitry Karpeev ierr = PetscSFSetGraph(bmsf,m,2*bm,PETSC_NULL,PETSC_COPY_VALUES,bmedges,PETSC_OWN_POINTER);CHKERRQ(ierr); 1456629c3df2SDmitry Karpeev ierr = PetscSFReduceBegin(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr); 1457629c3df2SDmitry Karpeev ierr = PetscSFReduceEnd(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr); 1458629c3df2SDmitry Karpeev ierr = PetscFree(bmdnnz);CHKERRQ(ierr); 1459629c3df2SDmitry Karpeev ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr); 1460629c3df2SDmitry Karpeev } 146122d28d08SBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1462629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1463629c3df2SDmitry Karpeev } 1464629c3df2SDmitry Karpeev ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr); 1465629c3df2SDmitry Karpeev ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr); 1466629c3df2SDmitry Karpeev ierr = PetscFree(dnnz);CHKERRQ(ierr); 1467629c3df2SDmitry Karpeev 1468629c3df2SDmitry Karpeev /* Fill by row */ 1469629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1470629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1471629c3df2SDmitry Karpeev IS bNis; 1472629c3df2SDmitry Karpeev PetscInt bN; 1473629c3df2SDmitry Karpeev const PetscInt *bNindices; 1474629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1475629c3df2SDmitry Karpeev ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr); 1476629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1477629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1478629c3df2SDmitry Karpeev Mat B; 1479629c3df2SDmitry Karpeev PetscInt bm, br; 1480629c3df2SDmitry Karpeev const PetscInt *bmindices; 1481629c3df2SDmitry Karpeev B = nest->m[i][j]; 1482629c3df2SDmitry Karpeev if (!B) continue; 1483629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1484629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1485629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1486629c3df2SDmitry Karpeev PetscInt row = bmindices[br], brncols, *cols; 1487629c3df2SDmitry Karpeev const PetscInt *brcols; 1488629c3df2SDmitry Karpeev const PetscScalar *brcoldata; 1489629c3df2SDmitry Karpeev ierr = MatGetRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1490629c3df2SDmitry Karpeev ierr = PetscMalloc(brncols*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1491629c3df2SDmitry Karpeev for (k=0; k<brncols; k++) 1492629c3df2SDmitry Karpeev cols[k] = bNindices[brcols[k]]; 1493629c3df2SDmitry Karpeev /* 1494629c3df2SDmitry Karpeev Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match. 1495629c3df2SDmitry Karpeev Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES. 1496629c3df2SDmitry Karpeev */ 1497a2ea699eSBarry Smith ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr); 1498629c3df2SDmitry Karpeev ierr = MatRestoreRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1499629c3df2SDmitry Karpeev ierr = PetscFree(cols);CHKERRQ(ierr); 1500629c3df2SDmitry Karpeev } 1501629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1502629c3df2SDmitry Karpeev } 1503a2ea699eSBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1504629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1505629c3df2SDmitry Karpeev } 1506629c3df2SDmitry Karpeev ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1507629c3df2SDmitry Karpeev ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1508629c3df2SDmitry Karpeev PetscFunctionReturn(0); 1509629c3df2SDmitry Karpeev } 1510629c3df2SDmitry Karpeev EXTERN_C_END 1511629c3df2SDmitry Karpeev 1512659c6bb0SJed Brown /*MC 1513659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1514659c6bb0SJed Brown 1515659c6bb0SJed Brown Level: intermediate 1516659c6bb0SJed Brown 1517659c6bb0SJed Brown Notes: 1518659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1519659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1520950540a4SJed Brown It is usually used with DMComposite and DMCreateMatrix() 1521659c6bb0SJed Brown 1522659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1523659c6bb0SJed Brown M*/ 1524c8883902SJed Brown EXTERN_C_BEGIN 1525c8883902SJed Brown #undef __FUNCT__ 1526c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 1527c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A) 1528c8883902SJed Brown { 1529c8883902SJed Brown Mat_Nest *s; 1530c8883902SJed Brown PetscErrorCode ierr; 1531c8883902SJed Brown 1532c8883902SJed Brown PetscFunctionBegin; 1533c8883902SJed Brown ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr); 1534c8883902SJed Brown A->data = (void*)s; 1535e7c19651SJed Brown 1536e7c19651SJed Brown s->nr = -1; 1537e7c19651SJed Brown s->nc = -1; 1538c8883902SJed Brown s->m = PETSC_NULL; 1539e7c19651SJed Brown s->splitassembly = PETSC_FALSE; 1540c8883902SJed Brown 1541c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 1542c8883902SJed Brown A->ops->mult = MatMult_Nest; 15439194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 1544c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 15459194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1546c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1547c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1548c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1549c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1550c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1551c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1552c8883902SJed Brown A->ops->view = MatView_Nest; 1553c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1554c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1555c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1556429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 1557429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 1558a061e289SJed Brown A->ops->scale = MatScale_Nest; 1559a061e289SJed Brown A->ops->shift = MatShift_Nest; 1560c8883902SJed Brown 1561c8883902SJed Brown A->spptr = 0; 1562c8883902SJed Brown A->same_nonzero = PETSC_FALSE; 1563c8883902SJed Brown A->assembled = PETSC_FALSE; 1564c8883902SJed Brown 1565c8883902SJed Brown /* expose Nest api's */ 1566c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C", "MatNestGetSubMat_Nest", MatNestGetSubMat_Nest);CHKERRQ(ierr); 15670782ca92SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C", "MatNestSetSubMat_Nest", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1568c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C", "MatNestGetSubMats_Nest", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1569c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C", "MatNestGetSize_Nest", MatNestGetSize_Nest);CHKERRQ(ierr); 1570900e7ff2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetISs_C", "MatNestGetISs_Nest", MatNestGetISs_Nest);CHKERRQ(ierr); 1571900e7ff2SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetLocalISs_C", "MatNestGetLocalISs_Nest", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 1572c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C", "MatNestSetVecType_Nest", MatNestSetVecType_Nest);CHKERRQ(ierr); 1573c8883902SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C", "MatNestSetSubMats_Nest", MatNestSetSubMats_Nest);CHKERRQ(ierr); 1574c8883902SJed Brown 1575c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1576c8883902SJed Brown PetscFunctionReturn(0); 1577c8883902SJed Brown } 157886e80854SHong Zhang EXTERN_C_END 1579