1d8588912SDave May 2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I "petscmat.h" I*/ 30c312b8eSJed Brown #include <petscsf.h> 4d8588912SDave May 5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]); 62a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left); 7c8883902SJed Brown 8d8588912SDave May /* private functions */ 9d8588912SDave May #undef __FUNCT__ 108188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private" 118188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N) 12d8588912SDave May { 13d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 148188e55aSJed Brown PetscInt i,j; 15d8588912SDave May PetscErrorCode ierr; 16d8588912SDave May 17d8588912SDave May PetscFunctionBegin; 188188e55aSJed Brown *m = *n = *M = *N = 0; 198188e55aSJed Brown for (i=0; i<bA->nr; i++) { /* rows */ 208188e55aSJed Brown PetscInt sm,sM; 218188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr); 228188e55aSJed Brown ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr); 238188e55aSJed Brown *m += sm; 248188e55aSJed Brown *M += sM; 25d8588912SDave May } 268188e55aSJed Brown for (j=0; j<bA->nc; j++) { /* cols */ 278188e55aSJed Brown PetscInt sn,sN; 288188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr); 298188e55aSJed Brown ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr); 308188e55aSJed Brown *n += sn; 318188e55aSJed Brown *N += sN; 32d8588912SDave May } 33d8588912SDave May PetscFunctionReturn(0); 34d8588912SDave May } 35d8588912SDave May 36d8588912SDave May /* operations */ 37d8588912SDave May #undef __FUNCT__ 38d8588912SDave May #define __FUNCT__ "MatMult_Nest" 39207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 40d8588912SDave May { 41d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 42207556f9SJed Brown Vec *bx = bA->right,*by = bA->left; 43207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 44d8588912SDave May PetscErrorCode ierr; 45d8588912SDave May 46d8588912SDave May PetscFunctionBegin; 47207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 48207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 49207556f9SJed Brown for (i=0; i<nr; i++) { 50d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 51207556f9SJed Brown for (j=0; j<nc; j++) { 52207556f9SJed Brown if (!bA->m[i][j]) continue; 53d8588912SDave May /* y[i] <- y[i] + A[i][j] * x[j] */ 54d8588912SDave May ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 55d8588912SDave May } 56d8588912SDave May } 57207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 58207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 59d8588912SDave May PetscFunctionReturn(0); 60d8588912SDave May } 61d8588912SDave May 62d8588912SDave May #undef __FUNCT__ 639194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest" 649194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z) 659194d70fSJed Brown { 669194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 679194d70fSJed Brown Vec *bx = bA->right,*bz = bA->left; 689194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 699194d70fSJed Brown PetscErrorCode ierr; 709194d70fSJed Brown 719194d70fSJed Brown PetscFunctionBegin; 729194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 739194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 749194d70fSJed Brown for (i=0; i<nr; i++) { 759194d70fSJed Brown if (y != z) { 769194d70fSJed Brown Vec by; 779194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 789194d70fSJed Brown ierr = VecCopy(by,bz[i]);CHKERRQ(ierr); 79336d21e7SJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 809194d70fSJed Brown } 819194d70fSJed Brown for (j=0; j<nc; j++) { 829194d70fSJed Brown if (!bA->m[i][j]) continue; 839194d70fSJed Brown /* y[i] <- y[i] + A[i][j] * x[j] */ 849194d70fSJed Brown ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr); 859194d70fSJed Brown } 869194d70fSJed Brown } 879194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 889194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 899194d70fSJed Brown PetscFunctionReturn(0); 909194d70fSJed Brown } 919194d70fSJed Brown 929194d70fSJed Brown #undef __FUNCT__ 93d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest" 94207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 95d8588912SDave May { 96d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 97207556f9SJed Brown Vec *bx = bA->left,*by = bA->right; 98207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 99d8588912SDave May PetscErrorCode ierr; 100d8588912SDave May 101d8588912SDave May PetscFunctionBegin; 102609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 103609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 104207556f9SJed Brown for (j=0; j<nc; j++) { 105609e31cbSJed Brown ierr = VecZeroEntries(by[j]);CHKERRQ(ierr); 106609e31cbSJed Brown for (i=0; i<nr; i++) { 1076c75ac25SJed Brown if (!bA->m[i][j]) continue; 108609e31cbSJed Brown /* y[j] <- y[j] + (A[i][j])^T * x[i] */ 109609e31cbSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr); 110d8588912SDave May } 111d8588912SDave May } 112609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 113609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 114d8588912SDave May PetscFunctionReturn(0); 115d8588912SDave May } 116d8588912SDave May 117d8588912SDave May #undef __FUNCT__ 1189194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest" 1199194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z) 1209194d70fSJed Brown { 1219194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 1229194d70fSJed Brown Vec *bx = bA->left,*bz = bA->right; 1239194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 1249194d70fSJed Brown PetscErrorCode ierr; 1259194d70fSJed Brown 1269194d70fSJed Brown PetscFunctionBegin; 1279194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 1289194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 1299194d70fSJed Brown for (j=0; j<nc; j++) { 1309194d70fSJed Brown if (y != z) { 1319194d70fSJed Brown Vec by; 1329194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 1339194d70fSJed Brown ierr = VecCopy(by,bz[j]);CHKERRQ(ierr); 1349194d70fSJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 1359194d70fSJed Brown } 1369194d70fSJed Brown for (i=0; i<nr; i++) { 1376c75ac25SJed Brown if (!bA->m[i][j]) continue; 1389194d70fSJed Brown /* z[j] <- y[j] + (A[i][j])^T * x[i] */ 1399194d70fSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr); 1409194d70fSJed Brown } 1419194d70fSJed Brown } 1429194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 1439194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 1449194d70fSJed Brown PetscFunctionReturn(0); 1459194d70fSJed Brown } 1469194d70fSJed Brown 1479194d70fSJed Brown #undef __FUNCT__ 148f8170845SAlex Fikl #define __FUNCT__ "MatTranspose_Nest" 149f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B) 150f8170845SAlex Fikl { 151f8170845SAlex Fikl Mat_Nest *bA = (Mat_Nest*)A->data, *bC; 152f8170845SAlex Fikl Mat C; 153f8170845SAlex Fikl PetscInt i,j,nr = bA->nr,nc = bA->nc; 154f8170845SAlex Fikl PetscErrorCode ierr; 155f8170845SAlex Fikl 156f8170845SAlex Fikl PetscFunctionBegin; 157f8170845SAlex Fikl if (reuse == MAT_REUSE_MATRIX && A == *B && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place"); 158f8170845SAlex Fikl 159f8170845SAlex Fikl if (reuse == MAT_INITIAL_MATRIX || A == *B) { 160f8170845SAlex Fikl Mat *subs; 161f8170845SAlex Fikl IS *is_row,*is_col; 162f8170845SAlex Fikl 163f8170845SAlex Fikl ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr); 164f8170845SAlex Fikl ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr); 165f8170845SAlex Fikl ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr); 166*ddeb9bd8SAlex Fikl if (reuse == MAT_REUSE_MATRIX) { 167*ddeb9bd8SAlex Fikl for (i=0; i<nr; i++) { 168*ddeb9bd8SAlex Fikl for (j=0; j<nc; j++) { 169*ddeb9bd8SAlex Fikl subs[i + nr * j] = bA->m[i][j]; 170*ddeb9bd8SAlex Fikl } 171*ddeb9bd8SAlex Fikl } 172*ddeb9bd8SAlex Fikl } 173*ddeb9bd8SAlex Fikl 174f8170845SAlex Fikl ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr); 175f8170845SAlex Fikl ierr = PetscFree(subs);CHKERRQ(ierr); 176f8170845SAlex Fikl ierr = PetscFree(is_row);CHKERRQ(ierr); 177f8170845SAlex Fikl ierr = PetscFree(is_col);CHKERRQ(ierr); 178f8170845SAlex Fikl } else { 179f8170845SAlex Fikl C = *B; 180f8170845SAlex Fikl } 181f8170845SAlex Fikl 182f8170845SAlex Fikl bC = (Mat_Nest*)C->data; 183f8170845SAlex Fikl for (i=0; i<nr; i++) { 184f8170845SAlex Fikl for (j=0; j<nc; j++) { 185f8170845SAlex Fikl if (bA->m[i][j]) { 186f8170845SAlex Fikl ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr); 187f8170845SAlex Fikl } else { 188f8170845SAlex Fikl bC->m[j][i] = NULL; 189f8170845SAlex Fikl } 190f8170845SAlex Fikl } 191f8170845SAlex Fikl } 192f8170845SAlex Fikl 193f8170845SAlex Fikl if (reuse == MAT_INITIAL_MATRIX || A != *B) { 194f8170845SAlex Fikl *B = C; 195f8170845SAlex Fikl } else { 196f8170845SAlex Fikl ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr); 197f8170845SAlex Fikl } 198f8170845SAlex Fikl PetscFunctionReturn(0); 199f8170845SAlex Fikl } 200f8170845SAlex Fikl 201f8170845SAlex Fikl #undef __FUNCT__ 202e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList" 203e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 204e2d7f03fSJed Brown { 205e2d7f03fSJed Brown PetscErrorCode ierr; 206e2d7f03fSJed Brown IS *lst = *list; 207e2d7f03fSJed Brown PetscInt i; 208e2d7f03fSJed Brown 209e2d7f03fSJed Brown PetscFunctionBegin; 210e2d7f03fSJed Brown if (!lst) PetscFunctionReturn(0); 2116bf464f9SBarry Smith for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);} 212e2d7f03fSJed Brown ierr = PetscFree(lst);CHKERRQ(ierr); 2130298fd71SBarry Smith *list = NULL; 214e2d7f03fSJed Brown PetscFunctionReturn(0); 215e2d7f03fSJed Brown } 216e2d7f03fSJed Brown 217e2d7f03fSJed Brown #undef __FUNCT__ 218d8588912SDave May #define __FUNCT__ "MatDestroy_Nest" 219207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A) 220d8588912SDave May { 221d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 222d8588912SDave May PetscInt i,j; 223d8588912SDave May PetscErrorCode ierr; 224d8588912SDave May 225d8588912SDave May PetscFunctionBegin; 226d8588912SDave May /* release the matrices and the place holders */ 227e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 228e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 229e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 230e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 231d8588912SDave May 232d8588912SDave May ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 233d8588912SDave May ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 234d8588912SDave May 235207556f9SJed Brown ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 236207556f9SJed Brown 237d8588912SDave May /* release the matrices and the place holders */ 238d8588912SDave May if (vs->m) { 239d8588912SDave May for (i=0; i<vs->nr; i++) { 240d8588912SDave May for (j=0; j<vs->nc; j++) { 2416bf464f9SBarry Smith ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr); 242d8588912SDave May } 243d8588912SDave May ierr = PetscFree(vs->m[i]);CHKERRQ(ierr); 244d8588912SDave May } 245d8588912SDave May ierr = PetscFree(vs->m);CHKERRQ(ierr); 246d8588912SDave May } 247bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 248d8588912SDave May 249bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr); 250bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr); 251bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr); 252bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr); 253bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr); 254bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr); 255bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr); 256bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr); 257d8588912SDave May PetscFunctionReturn(0); 258d8588912SDave May } 259d8588912SDave May 260d8588912SDave May #undef __FUNCT__ 261d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest" 262207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 263d8588912SDave May { 264d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 265d8588912SDave May PetscInt i,j; 266d8588912SDave May PetscErrorCode ierr; 267d8588912SDave May 268d8588912SDave May PetscFunctionBegin; 269d8588912SDave May for (i=0; i<vs->nr; i++) { 270d8588912SDave May for (j=0; j<vs->nc; j++) { 271e7c19651SJed Brown if (vs->m[i][j]) { 272e7c19651SJed Brown ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); 273e7c19651SJed Brown if (!vs->splitassembly) { 274e7c19651SJed Brown /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested 275e7c19651SJed 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 276e7c19651SJed Brown * already performing an assembly, but the result would by more complicated and appears to offer less 277e7c19651SJed Brown * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an 278e7c19651SJed Brown * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives". 279e7c19651SJed Brown */ 280e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 281e7c19651SJed Brown } 282e7c19651SJed Brown } 283d8588912SDave May } 284d8588912SDave May } 285d8588912SDave May PetscFunctionReturn(0); 286d8588912SDave May } 287d8588912SDave May 288d8588912SDave May #undef __FUNCT__ 289d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest" 290207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 291d8588912SDave May { 292d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 293d8588912SDave May PetscInt i,j; 294d8588912SDave May PetscErrorCode ierr; 295d8588912SDave May 296d8588912SDave May PetscFunctionBegin; 297d8588912SDave May for (i=0; i<vs->nr; i++) { 298d8588912SDave May for (j=0; j<vs->nc; j++) { 299e7c19651SJed Brown if (vs->m[i][j]) { 300e7c19651SJed Brown if (vs->splitassembly) { 301e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 302e7c19651SJed Brown } 303e7c19651SJed Brown } 304d8588912SDave May } 305d8588912SDave May } 306d8588912SDave May PetscFunctionReturn(0); 307d8588912SDave May } 308d8588912SDave May 309d8588912SDave May #undef __FUNCT__ 310f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow" 311f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 312d8588912SDave May { 313207556f9SJed Brown PetscErrorCode ierr; 314f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 315f349c1fdSJed Brown PetscInt j; 316f349c1fdSJed Brown Mat sub; 317d8588912SDave May 318d8588912SDave May PetscFunctionBegin; 3190298fd71SBarry Smith sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */ 320f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 3214994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 322f349c1fdSJed Brown *B = sub; 323f349c1fdSJed Brown PetscFunctionReturn(0); 324d8588912SDave May } 325d8588912SDave May 326f349c1fdSJed Brown #undef __FUNCT__ 327f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol" 328f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 329f349c1fdSJed Brown { 330207556f9SJed Brown PetscErrorCode ierr; 331f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 332f349c1fdSJed Brown PetscInt i; 333f349c1fdSJed Brown Mat sub; 334f349c1fdSJed Brown 335f349c1fdSJed Brown PetscFunctionBegin; 3360298fd71SBarry Smith sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */ 337f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 3384994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 339f349c1fdSJed Brown *B = sub; 340f349c1fdSJed Brown PetscFunctionReturn(0); 341d8588912SDave May } 342d8588912SDave May 343f349c1fdSJed Brown #undef __FUNCT__ 344f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS" 345f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 346f349c1fdSJed Brown { 347f349c1fdSJed Brown PetscErrorCode ierr; 348f349c1fdSJed Brown PetscInt i; 349f349c1fdSJed Brown PetscBool flg; 350f349c1fdSJed Brown 351f349c1fdSJed Brown PetscFunctionBegin; 352f349c1fdSJed Brown PetscValidPointer(list,3); 353f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 354f349c1fdSJed Brown PetscValidIntPointer(found,5); 355f349c1fdSJed Brown *found = -1; 356f349c1fdSJed Brown for (i=0; i<n; i++) { 357207556f9SJed Brown if (!list[i]) continue; 358f349c1fdSJed Brown ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr); 359f349c1fdSJed Brown if (flg) { 360f349c1fdSJed Brown *found = i; 361f349c1fdSJed Brown PetscFunctionReturn(0); 362f349c1fdSJed Brown } 363f349c1fdSJed Brown } 364ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set"); 365f349c1fdSJed Brown PetscFunctionReturn(0); 366f349c1fdSJed Brown } 367f349c1fdSJed Brown 368f349c1fdSJed Brown #undef __FUNCT__ 3698188e55aSJed Brown #define __FUNCT__ "MatNestGetRow" 3708188e55aSJed Brown /* Get a block row as a new MatNest */ 3718188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 3728188e55aSJed Brown { 3738188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3748188e55aSJed Brown char keyname[256]; 3758188e55aSJed Brown PetscErrorCode ierr; 3768188e55aSJed Brown 3778188e55aSJed Brown PetscFunctionBegin; 3780298fd71SBarry Smith *B = NULL; 3798caf3d72SBarry Smith ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr); 3808188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 3818188e55aSJed Brown if (*B) PetscFunctionReturn(0); 3828188e55aSJed Brown 383ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 38426fbe8dcSKarl Rupp 3858188e55aSJed Brown (*B)->assembled = A->assembled; 38626fbe8dcSKarl Rupp 3878188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 3888188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 3898188e55aSJed Brown PetscFunctionReturn(0); 3908188e55aSJed Brown } 3918188e55aSJed Brown 3928188e55aSJed Brown #undef __FUNCT__ 393f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat" 394f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 395f349c1fdSJed Brown { 396f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 3978188e55aSJed Brown PetscErrorCode ierr; 3986b3a5b13SJed Brown PetscInt row,col; 399e072481dSJed Brown PetscBool same,isFullCol,isFullColGlobal; 400f349c1fdSJed Brown 401f349c1fdSJed Brown PetscFunctionBegin; 4028188e55aSJed Brown /* Check if full column space. This is a hack */ 4038188e55aSJed Brown isFullCol = PETSC_FALSE; 404251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 4058188e55aSJed Brown if (same) { 40677019fcaSJed Brown PetscInt n,first,step,i,an,am,afirst,astep; 4078188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 4088188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 40977019fcaSJed Brown isFullCol = PETSC_TRUE; 41005ce4453SJed Brown for (i=0,an=A->cmap->rstart; i<vs->nc; i++) { 41177019fcaSJed Brown ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr); 41277019fcaSJed Brown ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr); 41377019fcaSJed Brown if (afirst != an || astep != step) isFullCol = PETSC_FALSE; 41477019fcaSJed Brown an += am; 41577019fcaSJed Brown } 41605ce4453SJed Brown if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE; 4178188e55aSJed Brown } 418b2566f29SBarry Smith ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr); 4198188e55aSJed Brown 420e072481dSJed Brown if (isFullColGlobal) { 4218188e55aSJed Brown PetscInt row; 4228188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 4238188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 4248188e55aSJed Brown } else { 425f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 426f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 427b6480e04SStefano Zampini if (!vs->m[row][col]) { 428b6480e04SStefano Zampini PetscInt lr,lc; 429b6480e04SStefano Zampini 430b6480e04SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr); 431b6480e04SStefano Zampini ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr); 432b6480e04SStefano Zampini ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr); 433b6480e04SStefano Zampini ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 434b6480e04SStefano Zampini ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr); 435b6480e04SStefano Zampini ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 436b6480e04SStefano Zampini ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 437b6480e04SStefano Zampini } 438f349c1fdSJed Brown *B = vs->m[row][col]; 4398188e55aSJed Brown } 440f349c1fdSJed Brown PetscFunctionReturn(0); 441f349c1fdSJed Brown } 442f349c1fdSJed Brown 443f349c1fdSJed Brown #undef __FUNCT__ 444f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest" 445f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 446f349c1fdSJed Brown { 447f349c1fdSJed Brown PetscErrorCode ierr; 448f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 449f349c1fdSJed Brown Mat sub; 450f349c1fdSJed Brown 451f349c1fdSJed Brown PetscFunctionBegin; 452f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 453f349c1fdSJed Brown switch (reuse) { 454f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 4557874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 456f349c1fdSJed Brown *B = sub; 457f349c1fdSJed Brown break; 458f349c1fdSJed Brown case MAT_REUSE_MATRIX: 459ce94432eSBarry Smith if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 460f349c1fdSJed Brown break; 461f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 462f349c1fdSJed Brown break; 463511c6705SHong Zhang case MAT_INPLACE_MATRIX: /* Nothing to do */ 464511c6705SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet"); 465511c6705SHong Zhang break; 466f349c1fdSJed Brown } 467f349c1fdSJed Brown PetscFunctionReturn(0); 468f349c1fdSJed Brown } 469f349c1fdSJed Brown 470f349c1fdSJed Brown #undef __FUNCT__ 471f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest" 472f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 473f349c1fdSJed Brown { 474f349c1fdSJed Brown PetscErrorCode ierr; 475f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 476f349c1fdSJed Brown Mat sub; 477f349c1fdSJed Brown 478f349c1fdSJed Brown PetscFunctionBegin; 479f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 480f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 481f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 482f349c1fdSJed Brown *B = sub; 483d8588912SDave May PetscFunctionReturn(0); 484d8588912SDave May } 485d8588912SDave May 486d8588912SDave May #undef __FUNCT__ 487d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest" 488207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 489d8588912SDave May { 490d8588912SDave May PetscErrorCode ierr; 491f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 492f349c1fdSJed Brown Mat sub; 493d8588912SDave May 494d8588912SDave May PetscFunctionBegin; 495f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 496ce94432eSBarry Smith if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 497f349c1fdSJed Brown if (sub) { 498ce94432eSBarry Smith if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 4996bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 500d8588912SDave May } 501d8588912SDave May PetscFunctionReturn(0); 502d8588912SDave May } 503d8588912SDave May 504d8588912SDave May #undef __FUNCT__ 5057874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest" 5067874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 5077874fa86SDave May { 5087874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 5097874fa86SDave May PetscInt i; 5107874fa86SDave May PetscErrorCode ierr; 5117874fa86SDave May 5127874fa86SDave May PetscFunctionBegin; 5137874fa86SDave May for (i=0; i<bA->nr; i++) { 514429bac76SJed Brown Vec bv; 515429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 5167874fa86SDave May if (bA->m[i][i]) { 517429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 5187874fa86SDave May } else { 5195159a857SMatthew G. Knepley ierr = VecSet(bv,0.0);CHKERRQ(ierr); 5207874fa86SDave May } 521429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 5227874fa86SDave May } 5237874fa86SDave May PetscFunctionReturn(0); 5247874fa86SDave May } 5257874fa86SDave May 5267874fa86SDave May #undef __FUNCT__ 5277874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest" 5287874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 5297874fa86SDave May { 5307874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 531429bac76SJed Brown Vec bl,*br; 5327874fa86SDave May PetscInt i,j; 5337874fa86SDave May PetscErrorCode ierr; 5347874fa86SDave May 5357874fa86SDave May PetscFunctionBegin; 5363f800ebeSJed Brown ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr); 5372e6472ebSElliott Sales de Andrade if (r) { 538429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 5392e6472ebSElliott Sales de Andrade } 5402e6472ebSElliott Sales de Andrade bl = NULL; 5417874fa86SDave May for (i=0; i<bA->nr; i++) { 5422e6472ebSElliott Sales de Andrade if (l) { 543429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 5442e6472ebSElliott Sales de Andrade } 5457874fa86SDave May for (j=0; j<bA->nc; j++) { 5467874fa86SDave May if (bA->m[i][j]) { 547429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 5487874fa86SDave May } 5497874fa86SDave May } 5502e6472ebSElliott Sales de Andrade if (l) { 551a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 5527874fa86SDave May } 5532e6472ebSElliott Sales de Andrade } 5542e6472ebSElliott Sales de Andrade if (r) { 555429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 5562e6472ebSElliott Sales de Andrade } 557429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 5587874fa86SDave May PetscFunctionReturn(0); 5597874fa86SDave May } 5607874fa86SDave May 5617874fa86SDave May #undef __FUNCT__ 562a061e289SJed Brown #define __FUNCT__ "MatScale_Nest" 563a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 564a061e289SJed Brown { 565a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 566a061e289SJed Brown PetscInt i,j; 567a061e289SJed Brown PetscErrorCode ierr; 568a061e289SJed Brown 569a061e289SJed Brown PetscFunctionBegin; 570a061e289SJed Brown for (i=0; i<bA->nr; i++) { 571a061e289SJed Brown for (j=0; j<bA->nc; j++) { 572a061e289SJed Brown if (bA->m[i][j]) { 573a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 574a061e289SJed Brown } 575a061e289SJed Brown } 576a061e289SJed Brown } 577a061e289SJed Brown PetscFunctionReturn(0); 578a061e289SJed Brown } 579a061e289SJed Brown 580a061e289SJed Brown #undef __FUNCT__ 581a061e289SJed Brown #define __FUNCT__ "MatShift_Nest" 582a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 583a061e289SJed Brown { 584a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 585a061e289SJed Brown PetscInt i; 586a061e289SJed Brown PetscErrorCode ierr; 587a061e289SJed Brown 588a061e289SJed Brown PetscFunctionBegin; 589a061e289SJed Brown for (i=0; i<bA->nr; i++) { 590ce94432eSBarry Smith if (!bA->m[i][i]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i); 591a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 592a061e289SJed Brown } 593a061e289SJed Brown PetscFunctionReturn(0); 594a061e289SJed Brown } 595a061e289SJed Brown 596a061e289SJed Brown #undef __FUNCT__ 59713135bc6SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Nest" 59813135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is) 59913135bc6SAlex Fikl { 60013135bc6SAlex Fikl Mat_Nest *bA = (Mat_Nest*)A->data; 60113135bc6SAlex Fikl PetscInt i; 60213135bc6SAlex Fikl PetscErrorCode ierr; 60313135bc6SAlex Fikl 60413135bc6SAlex Fikl PetscFunctionBegin; 60513135bc6SAlex Fikl for (i=0; i<bA->nr; i++) { 60613135bc6SAlex Fikl Vec bv; 60713135bc6SAlex Fikl ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 60813135bc6SAlex Fikl if (bA->m[i][i]) { 60913135bc6SAlex Fikl ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr); 61013135bc6SAlex Fikl } 61113135bc6SAlex Fikl ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 61213135bc6SAlex Fikl } 61313135bc6SAlex Fikl PetscFunctionReturn(0); 61413135bc6SAlex Fikl } 61513135bc6SAlex Fikl 61613135bc6SAlex Fikl #undef __FUNCT__ 617f8170845SAlex Fikl #define __FUNCT__ "MatSetRandom_Nest" 618f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx) 619f8170845SAlex Fikl { 620f8170845SAlex Fikl Mat_Nest *bA = (Mat_Nest*)A->data; 621f8170845SAlex Fikl PetscInt i,j; 622f8170845SAlex Fikl PetscErrorCode ierr; 623f8170845SAlex Fikl 624f8170845SAlex Fikl PetscFunctionBegin; 625f8170845SAlex Fikl for (i=0; i<bA->nr; i++) { 626f8170845SAlex Fikl for (j=0; j<bA->nc; j++) { 627f8170845SAlex Fikl if (bA->m[i][j]) { 628f8170845SAlex Fikl ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr); 629f8170845SAlex Fikl } 630f8170845SAlex Fikl } 631f8170845SAlex Fikl } 632f8170845SAlex Fikl PetscFunctionReturn(0); 633f8170845SAlex Fikl } 634f8170845SAlex Fikl 635f8170845SAlex Fikl #undef __FUNCT__ 6362a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest" 6372a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left) 638d8588912SDave May { 639d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 640d8588912SDave May Vec *L,*R; 641d8588912SDave May MPI_Comm comm; 642d8588912SDave May PetscInt i,j; 643d8588912SDave May PetscErrorCode ierr; 644d8588912SDave May 645d8588912SDave May PetscFunctionBegin; 646ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 647d8588912SDave May if (right) { 648d8588912SDave May /* allocate R */ 649854ce69bSBarry Smith ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr); 650d8588912SDave May /* Create the right vectors */ 651d8588912SDave May for (j=0; j<bA->nc; j++) { 652d8588912SDave May for (i=0; i<bA->nr; i++) { 653d8588912SDave May if (bA->m[i][j]) { 6542a7a6963SBarry Smith ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr); 655d8588912SDave May break; 656d8588912SDave May } 657d8588912SDave May } 6586c4ed002SBarry Smith if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 659d8588912SDave May } 660f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 661d8588912SDave May /* hand back control to the nest vector */ 662d8588912SDave May for (j=0; j<bA->nc; j++) { 6636bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 664d8588912SDave May } 665d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 666d8588912SDave May } 667d8588912SDave May 668d8588912SDave May if (left) { 669d8588912SDave May /* allocate L */ 670854ce69bSBarry Smith ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr); 671d8588912SDave May /* Create the left vectors */ 672d8588912SDave May for (i=0; i<bA->nr; i++) { 673d8588912SDave May for (j=0; j<bA->nc; j++) { 674d8588912SDave May if (bA->m[i][j]) { 6752a7a6963SBarry Smith ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr); 676d8588912SDave May break; 677d8588912SDave May } 678d8588912SDave May } 6796c4ed002SBarry Smith if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 680d8588912SDave May } 681d8588912SDave May 682f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 683d8588912SDave May for (i=0; i<bA->nr; i++) { 6846bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 685d8588912SDave May } 686d8588912SDave May 687d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 688d8588912SDave May } 689d8588912SDave May PetscFunctionReturn(0); 690d8588912SDave May } 691d8588912SDave May 692d8588912SDave May #undef __FUNCT__ 693d8588912SDave May #define __FUNCT__ "MatView_Nest" 694207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 695d8588912SDave May { 696d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 697d8588912SDave May PetscBool isascii; 698d8588912SDave May PetscInt i,j; 699d8588912SDave May PetscErrorCode ierr; 700d8588912SDave May 701d8588912SDave May PetscFunctionBegin; 702251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 703d8588912SDave May if (isascii) { 704d8588912SDave May 705d8588912SDave May PetscViewerASCIIPrintf(viewer,"Matrix object: \n"); 706d8588912SDave May PetscViewerASCIIPushTab(viewer); /* push0 */ 707d8588912SDave May PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc); 708d8588912SDave May 709d8588912SDave May PetscViewerASCIIPrintf(viewer,"MatNest structure: \n"); 710d8588912SDave May for (i=0; i<bA->nr; i++) { 711d8588912SDave May for (j=0; j<bA->nc; j++) { 71219fd82e9SBarry Smith MatType type; 713270f95d7SJed Brown char name[256] = "",prefix[256] = ""; 714d8588912SDave May PetscInt NR,NC; 715d8588912SDave May PetscBool isNest = PETSC_FALSE; 716d8588912SDave May 717d8588912SDave May if (!bA->m[i][j]) { 7180298fd71SBarry Smith PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j); 719d8588912SDave May continue; 720d8588912SDave May } 721d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 722d8588912SDave May ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr); 7238caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 7248caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 725251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 726d8588912SDave May 727270f95d7SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 728d8588912SDave May 729d8588912SDave May if (isNest) { 730270f95d7SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 731d8588912SDave May ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 732270f95d7SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 733d8588912SDave May } 734d8588912SDave May } 735d8588912SDave May } 736d8588912SDave May PetscViewerASCIIPopTab(viewer); /* pop0 */ 737d8588912SDave May } 738d8588912SDave May PetscFunctionReturn(0); 739d8588912SDave May } 740d8588912SDave May 741d8588912SDave May #undef __FUNCT__ 742d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest" 743207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 744d8588912SDave May { 745d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 746d8588912SDave May PetscInt i,j; 747d8588912SDave May PetscErrorCode ierr; 748d8588912SDave May 749d8588912SDave May PetscFunctionBegin; 750d8588912SDave May for (i=0; i<bA->nr; i++) { 751d8588912SDave May for (j=0; j<bA->nc; j++) { 752d8588912SDave May if (!bA->m[i][j]) continue; 753d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 754d8588912SDave May } 755d8588912SDave May } 756d8588912SDave May PetscFunctionReturn(0); 757d8588912SDave May } 758d8588912SDave May 759d8588912SDave May #undef __FUNCT__ 760c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest" 761c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str) 762c222c20dSDavid Ham { 763c222c20dSDavid Ham Mat_Nest *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data; 764c222c20dSDavid Ham PetscInt i,j,nr = bA->nr,nc = bA->nc; 765c222c20dSDavid Ham PetscErrorCode ierr; 766c222c20dSDavid Ham 767c222c20dSDavid Ham PetscFunctionBegin; 768c222c20dSDavid Ham if (nr != bB->nr || nc != bB->nc) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot copy a Mat_Nest of block size (%D,%D) to a Mat_Nest of block size (%D,%D)",bB->nr,bB->nc,nr,nc); 769c222c20dSDavid Ham for (i=0; i<nr; i++) { 770c222c20dSDavid Ham for (j=0; j<nc; j++) { 77146a2b97cSJed Brown if (bA->m[i][j] && bB->m[i][j]) { 772c222c20dSDavid Ham ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr); 77346a2b97cSJed Brown } else if (bA->m[i][j] || bB->m[i][j]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j); 774c222c20dSDavid Ham } 775c222c20dSDavid Ham } 776c222c20dSDavid Ham PetscFunctionReturn(0); 777c222c20dSDavid Ham } 778c222c20dSDavid Ham 779c222c20dSDavid Ham #undef __FUNCT__ 780d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest" 781207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 782d8588912SDave May { 783d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 784841e96a3SJed Brown Mat *b; 785841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 786d8588912SDave May PetscErrorCode ierr; 787d8588912SDave May 788d8588912SDave May PetscFunctionBegin; 789785e854fSJed Brown ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr); 790841e96a3SJed Brown for (i=0; i<nr; i++) { 791841e96a3SJed Brown for (j=0; j<nc; j++) { 792841e96a3SJed Brown if (bA->m[i][j]) { 793841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 794841e96a3SJed Brown } else { 7950298fd71SBarry Smith b[i*nc+j] = NULL; 796d8588912SDave May } 797d8588912SDave May } 798d8588912SDave May } 799ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 800841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 801841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 8026bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 803d8588912SDave May } 804d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 805d8588912SDave May 806841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 807841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 808d8588912SDave May PetscFunctionReturn(0); 809d8588912SDave May } 810d8588912SDave May 811d8588912SDave May /* nest api */ 812d8588912SDave May #undef __FUNCT__ 813d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest" 814d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 815d8588912SDave May { 816d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 8175fd66863SKarl Rupp 818d8588912SDave May PetscFunctionBegin; 819ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 820ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 821d8588912SDave May *mat = bA->m[idxm][jdxm]; 822d8588912SDave May PetscFunctionReturn(0); 823d8588912SDave May } 824d8588912SDave May 825d8588912SDave May #undef __FUNCT__ 826d8588912SDave May #define __FUNCT__ "MatNestGetSubMat" 8279ba0d327SJed Brown /*@ 828d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 829d8588912SDave May 830d8588912SDave May Not collective 831d8588912SDave May 832d8588912SDave May Input Parameters: 833629881c0SJed Brown + A - nest matrix 834d8588912SDave May . idxm - index of the matrix within the nest matrix 835629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 836d8588912SDave May 837d8588912SDave May Output Parameter: 838d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 839d8588912SDave May 840d8588912SDave May Level: developer 841d8588912SDave May 842d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats() 843d8588912SDave May @*/ 8447087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 845d8588912SDave May { 846699a902aSJed Brown PetscErrorCode ierr; 847d8588912SDave May 848d8588912SDave May PetscFunctionBegin; 849699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 850d8588912SDave May PetscFunctionReturn(0); 851d8588912SDave May } 852d8588912SDave May 853d8588912SDave May #undef __FUNCT__ 8540782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest" 8550782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 8560782ca92SJed Brown { 8570782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 8580782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 8590782ca92SJed Brown PetscErrorCode ierr; 8600782ca92SJed Brown 8610782ca92SJed Brown PetscFunctionBegin; 862ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 863ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 8640782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 8650782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 8660782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 8670782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 8680782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 8690782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 870ce94432eSBarry Smith if (M != Mi || N != Ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni); 871ce94432eSBarry Smith if (m != mi || n != ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni); 87226fbe8dcSKarl Rupp 8730782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 8740782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 8750782ca92SJed Brown bA->m[idxm][jdxm] = mat; 8760782ca92SJed Brown PetscFunctionReturn(0); 8770782ca92SJed Brown } 8780782ca92SJed Brown 8790782ca92SJed Brown #undef __FUNCT__ 8800782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat" 8819ba0d327SJed Brown /*@ 8820782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 8830782ca92SJed Brown 8840782ca92SJed Brown Logically collective on the submatrix communicator 8850782ca92SJed Brown 8860782ca92SJed Brown Input Parameters: 8870782ca92SJed Brown + A - nest matrix 8880782ca92SJed Brown . idxm - index of the matrix within the nest matrix 8890782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 8900782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 8910782ca92SJed Brown 8920782ca92SJed Brown Notes: 8930782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 8940782ca92SJed Brown 8950782ca92SJed Brown This increments the reference count of the submatrix. 8960782ca92SJed Brown 8970782ca92SJed Brown Level: developer 8980782ca92SJed Brown 8990782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat() 9000782ca92SJed Brown @*/ 9010782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 9020782ca92SJed Brown { 9030782ca92SJed Brown PetscErrorCode ierr; 9040782ca92SJed Brown 9050782ca92SJed Brown PetscFunctionBegin; 9060782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 9070782ca92SJed Brown PetscFunctionReturn(0); 9080782ca92SJed Brown } 9090782ca92SJed Brown 9100782ca92SJed Brown #undef __FUNCT__ 911d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest" 912d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 913d8588912SDave May { 914d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 9155fd66863SKarl Rupp 916d8588912SDave May PetscFunctionBegin; 91726fbe8dcSKarl Rupp if (M) *M = bA->nr; 91826fbe8dcSKarl Rupp if (N) *N = bA->nc; 91926fbe8dcSKarl Rupp if (mat) *mat = bA->m; 920d8588912SDave May PetscFunctionReturn(0); 921d8588912SDave May } 922d8588912SDave May 923d8588912SDave May #undef __FUNCT__ 924d8588912SDave May #define __FUNCT__ "MatNestGetSubMats" 925d8588912SDave May /*@C 926d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 927d8588912SDave May 928d8588912SDave May Not collective 929d8588912SDave May 930d8588912SDave May Input Parameters: 931629881c0SJed Brown . A - nest matrix 932d8588912SDave May 933d8588912SDave May Output Parameter: 934629881c0SJed Brown + M - number of rows in the nest matrix 935d8588912SDave May . N - number of cols in the nest matrix 936629881c0SJed Brown - mat - 2d array of matrices 937d8588912SDave May 938d8588912SDave May Notes: 939d8588912SDave May 940d8588912SDave May The user should not free the array mat. 941d8588912SDave May 942351962e3SVincent Le Chenadec In Fortran, this routine has a calling sequence 943351962e3SVincent Le Chenadec $ call MatNestGetSubMats(A, M, N, mat, ierr) 944351962e3SVincent Le Chenadec where the space allocated for the optional argument mat is assumed large enough (if provided). 945351962e3SVincent Le Chenadec 946d8588912SDave May Level: developer 947d8588912SDave May 948d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat() 949d8588912SDave May @*/ 9507087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 951d8588912SDave May { 952699a902aSJed Brown PetscErrorCode ierr; 953d8588912SDave May 954d8588912SDave May PetscFunctionBegin; 955699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 956d8588912SDave May PetscFunctionReturn(0); 957d8588912SDave May } 958d8588912SDave May 959d8588912SDave May #undef __FUNCT__ 960d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest" 9617087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 962d8588912SDave May { 963d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 964d8588912SDave May 965d8588912SDave May PetscFunctionBegin; 96626fbe8dcSKarl Rupp if (M) *M = bA->nr; 96726fbe8dcSKarl Rupp if (N) *N = bA->nc; 968d8588912SDave May PetscFunctionReturn(0); 969d8588912SDave May } 970d8588912SDave May 971d8588912SDave May #undef __FUNCT__ 972d8588912SDave May #define __FUNCT__ "MatNestGetSize" 9739ba0d327SJed Brown /*@ 974d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 975d8588912SDave May 976d8588912SDave May Not collective 977d8588912SDave May 978d8588912SDave May Input Parameters: 979d8588912SDave May . A - nest matrix 980d8588912SDave May 981d8588912SDave May Output Parameter: 982629881c0SJed Brown + M - number of rows in the nested mat 983629881c0SJed Brown - N - number of cols in the nested mat 984d8588912SDave May 985d8588912SDave May Notes: 986d8588912SDave May 987d8588912SDave May Level: developer 988d8588912SDave May 989d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats() 990d8588912SDave May @*/ 9917087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 992d8588912SDave May { 993699a902aSJed Brown PetscErrorCode ierr; 994d8588912SDave May 995d8588912SDave May PetscFunctionBegin; 996699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 997d8588912SDave May PetscFunctionReturn(0); 998d8588912SDave May } 999d8588912SDave May 1000900e7ff2SJed Brown #undef __FUNCT__ 1001900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest" 1002f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 1003900e7ff2SJed Brown { 1004900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 1005900e7ff2SJed Brown PetscInt i; 1006900e7ff2SJed Brown 1007900e7ff2SJed Brown PetscFunctionBegin; 1008900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 1009900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 1010900e7ff2SJed Brown PetscFunctionReturn(0); 1011900e7ff2SJed Brown } 1012900e7ff2SJed Brown 1013900e7ff2SJed Brown #undef __FUNCT__ 1014900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs" 10153a4d7b9aSSatish Balay /*@C 1016900e7ff2SJed Brown MatNestGetISs - Returns the index sets partitioning the row and column spaces 1017900e7ff2SJed Brown 1018900e7ff2SJed Brown Not collective 1019900e7ff2SJed Brown 1020900e7ff2SJed Brown Input Parameters: 1021900e7ff2SJed Brown . A - nest matrix 1022900e7ff2SJed Brown 1023900e7ff2SJed Brown Output Parameter: 1024900e7ff2SJed Brown + rows - array of row index sets 1025900e7ff2SJed Brown - cols - array of column index sets 1026900e7ff2SJed Brown 1027900e7ff2SJed Brown Level: advanced 1028900e7ff2SJed Brown 1029900e7ff2SJed Brown Notes: 1030900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 1031900e7ff2SJed Brown 1032900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs() 1033900e7ff2SJed Brown @*/ 1034900e7ff2SJed Brown PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 1035900e7ff2SJed Brown { 1036900e7ff2SJed Brown PetscErrorCode ierr; 1037900e7ff2SJed Brown 1038900e7ff2SJed Brown PetscFunctionBegin; 1039900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1040900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 1041900e7ff2SJed Brown PetscFunctionReturn(0); 1042900e7ff2SJed Brown } 1043900e7ff2SJed Brown 1044900e7ff2SJed Brown #undef __FUNCT__ 1045900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest" 1046f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 1047900e7ff2SJed Brown { 1048900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 1049900e7ff2SJed Brown PetscInt i; 1050900e7ff2SJed Brown 1051900e7ff2SJed Brown PetscFunctionBegin; 1052900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 1053900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 1054900e7ff2SJed Brown PetscFunctionReturn(0); 1055900e7ff2SJed Brown } 1056900e7ff2SJed Brown 1057900e7ff2SJed Brown #undef __FUNCT__ 1058900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs" 1059900e7ff2SJed Brown /*@C 1060900e7ff2SJed Brown MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 1061900e7ff2SJed Brown 1062900e7ff2SJed Brown Not collective 1063900e7ff2SJed Brown 1064900e7ff2SJed Brown Input Parameters: 1065900e7ff2SJed Brown . A - nest matrix 1066900e7ff2SJed Brown 1067900e7ff2SJed Brown Output Parameter: 10680298fd71SBarry Smith + rows - array of row index sets (or NULL to ignore) 10690298fd71SBarry Smith - cols - array of column index sets (or NULL to ignore) 1070900e7ff2SJed Brown 1071900e7ff2SJed Brown Level: advanced 1072900e7ff2SJed Brown 1073900e7ff2SJed Brown Notes: 1074900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 1075900e7ff2SJed Brown 1076900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs() 1077900e7ff2SJed Brown @*/ 1078900e7ff2SJed Brown PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 1079900e7ff2SJed Brown { 1080900e7ff2SJed Brown PetscErrorCode ierr; 1081900e7ff2SJed Brown 1082900e7ff2SJed Brown PetscFunctionBegin; 1083900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1084900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 1085900e7ff2SJed Brown PetscFunctionReturn(0); 1086900e7ff2SJed Brown } 1087900e7ff2SJed Brown 1088207556f9SJed Brown #undef __FUNCT__ 1089207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest" 109019fd82e9SBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,VecType vtype) 1091207556f9SJed Brown { 1092207556f9SJed Brown PetscErrorCode ierr; 1093207556f9SJed Brown PetscBool flg; 1094207556f9SJed Brown 1095207556f9SJed Brown PetscFunctionBegin; 1096207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 1097207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 10982a7a6963SBarry Smith if (flg) A->ops->getvecs = MatCreateVecs_Nest; 109912b53f24SSatish Balay else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0; 1100207556f9SJed Brown PetscFunctionReturn(0); 1101207556f9SJed Brown } 1102207556f9SJed Brown 1103207556f9SJed Brown #undef __FUNCT__ 1104207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType" 1105207556f9SJed Brown /*@C 11062a7a6963SBarry Smith MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs() 1107207556f9SJed Brown 1108207556f9SJed Brown Not collective 1109207556f9SJed Brown 1110207556f9SJed Brown Input Parameters: 1111207556f9SJed Brown + A - nest matrix 1112207556f9SJed Brown - vtype - type to use for creating vectors 1113207556f9SJed Brown 1114207556f9SJed Brown Notes: 1115207556f9SJed Brown 1116207556f9SJed Brown Level: developer 1117207556f9SJed Brown 11182a7a6963SBarry Smith .seealso: MatCreateVecs() 1119207556f9SJed Brown @*/ 112019fd82e9SBarry Smith PetscErrorCode MatNestSetVecType(Mat A,VecType vtype) 1121207556f9SJed Brown { 1122207556f9SJed Brown PetscErrorCode ierr; 1123207556f9SJed Brown 1124207556f9SJed Brown PetscFunctionBegin; 112519fd82e9SBarry Smith ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr); 1126207556f9SJed Brown PetscFunctionReturn(0); 1127207556f9SJed Brown } 1128207556f9SJed Brown 1129d8588912SDave May #undef __FUNCT__ 1130c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest" 1131c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1132d8588912SDave May { 1133c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 1134c8883902SJed Brown PetscInt i,j,m,n,M,N; 1135d8588912SDave May PetscErrorCode ierr; 1136d8588912SDave May 1137d8588912SDave May PetscFunctionBegin; 1138c8883902SJed Brown s->nr = nr; 1139c8883902SJed Brown s->nc = nc; 1140d8588912SDave May 1141c8883902SJed Brown /* Create space for submatrices */ 1142854ce69bSBarry Smith ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr); 1143c8883902SJed Brown for (i=0; i<nr; i++) { 1144854ce69bSBarry Smith ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr); 1145d8588912SDave May } 1146c8883902SJed Brown for (i=0; i<nr; i++) { 1147c8883902SJed Brown for (j=0; j<nc; j++) { 1148c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1149c8883902SJed Brown if (a[i*nc+j]) { 1150c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1151d8588912SDave May } 1152d8588912SDave May } 1153d8588912SDave May } 1154d8588912SDave May 11558188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1156d8588912SDave May 1157854ce69bSBarry Smith ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr); 1158854ce69bSBarry Smith ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr); 1159c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1160c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1161d8588912SDave May 11628188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1163d8588912SDave May 1164c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1165c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1166c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1167c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1168c8883902SJed Brown 1169c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1170c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1171c8883902SJed Brown 11721795a4d1SJed Brown ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr); 1173d8588912SDave May PetscFunctionReturn(0); 1174d8588912SDave May } 1175d8588912SDave May 1176c8883902SJed Brown #undef __FUNCT__ 1177c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats" 1178c8883902SJed Brown /*@ 1179c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1180c8883902SJed Brown 1181c8883902SJed Brown Collective on Mat 1182c8883902SJed Brown 1183c8883902SJed Brown Input Parameter: 1184c8883902SJed Brown + N - nested matrix 1185c8883902SJed Brown . nr - number of nested row blocks 11860298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1187c8883902SJed Brown . nc - number of nested column blocks 11880298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 11890298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1190c8883902SJed Brown 1191c8883902SJed Brown Level: advanced 1192c8883902SJed Brown 1193c8883902SJed Brown .seealso: MatCreateNest(), MATNEST 1194c8883902SJed Brown @*/ 1195c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1196c8883902SJed Brown { 1197c8883902SJed Brown PetscErrorCode ierr; 1198c8883902SJed Brown PetscInt i; 1199c8883902SJed Brown 1200c8883902SJed Brown PetscFunctionBegin; 1201c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1202ce94432eSBarry Smith if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1203c8883902SJed Brown if (nr && is_row) { 1204c8883902SJed Brown PetscValidPointer(is_row,3); 1205c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1206c8883902SJed Brown } 1207ce94432eSBarry Smith if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 12081664e352SJed Brown if (nc && is_col) { 1209c8883902SJed Brown PetscValidPointer(is_col,5); 12109b30a8f6SBarry Smith for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1211c8883902SJed Brown } 1212c8883902SJed Brown if (nr*nc) PetscValidPointer(a,6); 1213c8883902SJed 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); 1214c8883902SJed Brown PetscFunctionReturn(0); 1215c8883902SJed Brown } 1216d8588912SDave May 121777019fcaSJed Brown #undef __FUNCT__ 121877019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private" 121945b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog) 122077019fcaSJed Brown { 122177019fcaSJed Brown PetscErrorCode ierr; 122277019fcaSJed Brown PetscBool flg; 122377019fcaSJed Brown PetscInt i,j,m,mi,*ix; 122477019fcaSJed Brown 122577019fcaSJed Brown PetscFunctionBegin; 122677019fcaSJed Brown for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 122777019fcaSJed Brown if (islocal[i]) { 122877019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 122977019fcaSJed Brown flg = PETSC_TRUE; /* We found a non-trivial entry */ 123077019fcaSJed Brown } else { 123177019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 123277019fcaSJed Brown } 123377019fcaSJed Brown m += mi; 123477019fcaSJed Brown } 123577019fcaSJed Brown if (flg) { 1236785e854fSJed Brown ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr); 123777019fcaSJed Brown for (i=0,n=0; i<n; i++) { 12380298fd71SBarry Smith ISLocalToGlobalMapping smap = NULL; 123977019fcaSJed Brown VecScatter scat; 124077019fcaSJed Brown IS isreq; 124177019fcaSJed Brown Vec lvec,gvec; 12423361c9a7SJed Brown union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x; 124377019fcaSJed Brown Mat sub; 124477019fcaSJed Brown 1245ce94432eSBarry Smith if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers"); 124677019fcaSJed Brown if (colflg) { 124777019fcaSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 124877019fcaSJed Brown } else { 124977019fcaSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 125077019fcaSJed Brown } 12510298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);} 125277019fcaSJed Brown if (islocal[i]) { 125377019fcaSJed Brown ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr); 125477019fcaSJed Brown } else { 125577019fcaSJed Brown ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr); 125677019fcaSJed Brown } 125777019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = j; 125877019fcaSJed Brown if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);} 125977019fcaSJed Brown /* 126077019fcaSJed Brown Now we need to extract the monolithic global indices that correspond to the given split global indices. 126177019fcaSJed 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. 126277019fcaSJed Brown The approach here is ugly because it uses VecScatter to move indices. 126377019fcaSJed Brown */ 126477019fcaSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr); 126577019fcaSJed Brown ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr); 126677019fcaSJed Brown ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr); 12670298fd71SBarry Smith ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr); 126877019fcaSJed Brown ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 126977019fcaSJed Brown for (j=0; j<mi; j++) x[j].integer = ix[m+j]; 127077019fcaSJed Brown ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr); 127177019fcaSJed Brown ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 127277019fcaSJed Brown ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 127377019fcaSJed Brown ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 127477019fcaSJed Brown for (j=0; j<mi; j++) ix[m+j] = x[j].integer; 127577019fcaSJed Brown ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr); 127677019fcaSJed Brown ierr = VecDestroy(&lvec);CHKERRQ(ierr); 127777019fcaSJed Brown ierr = VecDestroy(&gvec);CHKERRQ(ierr); 127877019fcaSJed Brown ierr = ISDestroy(&isreq);CHKERRQ(ierr); 127977019fcaSJed Brown ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 128077019fcaSJed Brown m += mi; 128177019fcaSJed Brown } 1282f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 128377019fcaSJed Brown } else { 12840298fd71SBarry Smith *ltog = NULL; 128577019fcaSJed Brown } 128677019fcaSJed Brown PetscFunctionReturn(0); 128777019fcaSJed Brown } 128877019fcaSJed Brown 128977019fcaSJed Brown 1290d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1291d8588912SDave May /* 1292d8588912SDave May nprocessors = NP 1293d8588912SDave May Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1)) 1294d8588912SDave May proc 0: => (g_0,h_0,) 1295d8588912SDave May proc 1: => (g_1,h_1,) 1296d8588912SDave May ... 1297d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1298d8588912SDave May 1299d8588912SDave May proc 0: proc 1: proc nprocs-1: 1300d8588912SDave May is[0] = (0,1,2,...,nlocal(g_0)-1) (0,1,...,nlocal(g_1)-1) (0,1,...,nlocal(g_NP-1)) 1301d8588912SDave May 1302d8588912SDave May proc 0: 1303d8588912SDave May is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1) 1304d8588912SDave May proc 1: 1305d8588912SDave May is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1) 1306d8588912SDave May 1307d8588912SDave May proc NP-1: 1308d8588912SDave May is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1) 1309d8588912SDave May */ 1310d8588912SDave May #undef __FUNCT__ 1311d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private" 1312841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1313d8588912SDave May { 1314e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 13158188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1316d8588912SDave May PetscErrorCode ierr; 13170298fd71SBarry Smith Mat sub = NULL; 1318d8588912SDave May 1319d8588912SDave May PetscFunctionBegin; 1320854ce69bSBarry Smith ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr); 1321854ce69bSBarry Smith ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr); 1322d8588912SDave May if (is_row) { /* valid IS is passed in */ 1323d8588912SDave May /* refs on is[] are incremeneted */ 1324e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1325d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 132626fbe8dcSKarl Rupp 1327e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1328d8588912SDave May } 13292ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 13308188e55aSJed Brown nsum = 0; 13318188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 13328188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1333ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 13340298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 1335ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 13368188e55aSJed Brown nsum += n; 13378188e55aSJed Brown } 1338ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 133930bc264bSJed Brown offset -= nsum; 1340e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1341f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 13420298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 13432ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1344ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1345e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 13462ae74bdbSJed Brown offset += n; 1347d8588912SDave May } 1348d8588912SDave May } 1349d8588912SDave May 1350d8588912SDave May if (is_col) { /* valid IS is passed in */ 1351d8588912SDave May /* refs on is[] are incremeneted */ 1352e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1353d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 135426fbe8dcSKarl Rupp 1355e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1356d8588912SDave May } 13572ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 13582ae74bdbSJed Brown offset = A->cmap->rstart; 13598188e55aSJed Brown nsum = 0; 13608188e55aSJed Brown for (j=0; j<vs->nc; j++) { 13618188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1362ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 13630298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 1364ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 13658188e55aSJed Brown nsum += n; 13668188e55aSJed Brown } 1367ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 136830bc264bSJed Brown offset -= nsum; 1369e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1370f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 13710298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 13722ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1373ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1374e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 13752ae74bdbSJed Brown offset += n; 1376d8588912SDave May } 1377d8588912SDave May } 1378e2d7f03fSJed Brown 1379e2d7f03fSJed Brown /* Set up the local ISs */ 1380785e854fSJed Brown ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 1381785e854fSJed Brown ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 1382e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1383e2d7f03fSJed Brown IS isloc; 13840298fd71SBarry Smith ISLocalToGlobalMapping rmap = NULL; 1385e2d7f03fSJed Brown PetscInt nlocal,bs; 1386e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 13870298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);} 1388207556f9SJed Brown if (rmap) { 1389e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1390e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1391e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1392e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1393207556f9SJed Brown } else { 1394207556f9SJed Brown nlocal = 0; 13950298fd71SBarry Smith isloc = NULL; 1396207556f9SJed Brown } 1397e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1398e2d7f03fSJed Brown offset += nlocal; 1399e2d7f03fSJed Brown } 14008188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1401e2d7f03fSJed Brown IS isloc; 14020298fd71SBarry Smith ISLocalToGlobalMapping cmap = NULL; 1403e2d7f03fSJed Brown PetscInt nlocal,bs; 1404e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 14050298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);} 1406207556f9SJed Brown if (cmap) { 1407e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1408e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1409e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1410e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1411207556f9SJed Brown } else { 1412207556f9SJed Brown nlocal = 0; 14130298fd71SBarry Smith isloc = NULL; 1414207556f9SJed Brown } 1415e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1416e2d7f03fSJed Brown offset += nlocal; 1417e2d7f03fSJed Brown } 14180189643fSJed Brown 141977019fcaSJed Brown /* Set up the aggregate ISLocalToGlobalMapping */ 142077019fcaSJed Brown { 142145b6f7e9SBarry Smith ISLocalToGlobalMapping rmap,cmap; 142245b6f7e9SBarry Smith ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr); 142345b6f7e9SBarry Smith ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr); 142477019fcaSJed Brown if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 142577019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 142677019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 142777019fcaSJed Brown } 142877019fcaSJed Brown 14290189643fSJed Brown #if defined(PETSC_USE_DEBUG) 14300189643fSJed Brown for (i=0; i<vs->nr; i++) { 14310189643fSJed Brown for (j=0; j<vs->nc; j++) { 14320189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 14330189643fSJed Brown Mat B = vs->m[i][j]; 14340189643fSJed Brown if (!B) continue; 14350189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 14360189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 14370189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 14380189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 14390189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 14400189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 1441ce94432eSBarry Smith if (M != Mi || N != Ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni); 1442ce94432eSBarry Smith if (m != mi || n != ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni); 14430189643fSJed Brown } 14440189643fSJed Brown } 14450189643fSJed Brown #endif 1446a061e289SJed Brown 1447a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1448a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1449a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1450a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1451a061e289SJed Brown } 1452a061e289SJed Brown } 1453a061e289SJed Brown A->assembled = PETSC_TRUE; 1454d8588912SDave May PetscFunctionReturn(0); 1455d8588912SDave May } 1456d8588912SDave May 1457d8588912SDave May #undef __FUNCT__ 1458d8588912SDave May #define __FUNCT__ "MatCreateNest" 145945c38901SJed Brown /*@C 1460659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1461659c6bb0SJed Brown 1462659c6bb0SJed Brown Collective on Mat 1463659c6bb0SJed Brown 1464659c6bb0SJed Brown Input Parameter: 1465659c6bb0SJed Brown + comm - Communicator for the new Mat 1466659c6bb0SJed Brown . nr - number of nested row blocks 14670298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1468659c6bb0SJed Brown . nc - number of nested column blocks 14690298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 14700298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1471659c6bb0SJed Brown 1472659c6bb0SJed Brown Output Parameter: 1473659c6bb0SJed Brown . B - new matrix 1474659c6bb0SJed Brown 1475659c6bb0SJed Brown Level: advanced 1476659c6bb0SJed Brown 1477950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST 1478659c6bb0SJed Brown @*/ 14797087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1480d8588912SDave May { 1481d8588912SDave May Mat A; 1482d8588912SDave May PetscErrorCode ierr; 1483d8588912SDave May 1484d8588912SDave May PetscFunctionBegin; 1485c8883902SJed Brown *B = 0; 1486d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1487c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 148891a28eb3SBarry Smith A->preallocated = PETSC_TRUE; 1489c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1490d8588912SDave May *B = A; 1491d8588912SDave May PetscFunctionReturn(0); 1492d8588912SDave May } 1493659c6bb0SJed Brown 1494629c3df2SDmitry Karpeev #undef __FUNCT__ 1495629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ" 1496cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1497629c3df2SDmitry Karpeev { 1498629c3df2SDmitry Karpeev PetscErrorCode ierr; 1499629c3df2SDmitry Karpeev Mat_Nest *nest = (Mat_Nest*)A->data; 150083b1a929SMark Adams PetscInt m,n,M,N,i,j,k,*dnnz,*onnz,rstart; 1501649b366bSFande Kong PetscInt cstart,cend; 1502629c3df2SDmitry Karpeev Mat C; 1503629c3df2SDmitry Karpeev 1504629c3df2SDmitry Karpeev PetscFunctionBegin; 1505629c3df2SDmitry Karpeev ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1506629c3df2SDmitry Karpeev ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1507649b366bSFande Kong ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); 1508629c3df2SDmitry Karpeev switch (reuse) { 1509629c3df2SDmitry Karpeev case MAT_INITIAL_MATRIX: 1510ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 1511629c3df2SDmitry Karpeev ierr = MatSetType(C,newtype);CHKERRQ(ierr); 1512629c3df2SDmitry Karpeev ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 1513629c3df2SDmitry Karpeev *newmat = C; 1514629c3df2SDmitry Karpeev break; 1515629c3df2SDmitry Karpeev case MAT_REUSE_MATRIX: 1516629c3df2SDmitry Karpeev C = *newmat; 1517629c3df2SDmitry Karpeev break; 1518ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse"); 1519629c3df2SDmitry Karpeev } 1520785e854fSJed Brown ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr); 1521629c3df2SDmitry Karpeev onnz = dnnz + m; 1522629c3df2SDmitry Karpeev for (k=0; k<m; k++) { 1523629c3df2SDmitry Karpeev dnnz[k] = 0; 1524629c3df2SDmitry Karpeev onnz[k] = 0; 1525629c3df2SDmitry Karpeev } 1526629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1527629c3df2SDmitry Karpeev IS bNis; 1528629c3df2SDmitry Karpeev PetscInt bN; 1529629c3df2SDmitry Karpeev const PetscInt *bNindices; 1530629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1531629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1532629c3df2SDmitry Karpeev ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr); 1533629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1534629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1535629c3df2SDmitry Karpeev PetscSF bmsf; 1536649b366bSFande Kong PetscSFNode *iremote; 1537629c3df2SDmitry Karpeev Mat B; 1538649b366bSFande Kong PetscInt bm, *sub_dnnz,*sub_onnz, br; 1539629c3df2SDmitry Karpeev const PetscInt *bmindices; 1540629c3df2SDmitry Karpeev B = nest->m[i][j]; 1541629c3df2SDmitry Karpeev if (!B) continue; 1542629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1543629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1544ce94432eSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr); 1545649b366bSFande Kong ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr); 1546649b366bSFande Kong ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr); 1547649b366bSFande Kong ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr); 1548649b366bSFande Kong for (k = 0; k < bm; ++k){ 1549649b366bSFande Kong sub_dnnz[k] = 0; 1550649b366bSFande Kong sub_onnz[k] = 0; 1551649b366bSFande Kong } 1552629c3df2SDmitry Karpeev /* 1553629c3df2SDmitry Karpeev Locate the owners for all of the locally-owned global row indices for this row block. 1554629c3df2SDmitry Karpeev These determine the roots of PetscSF used to communicate preallocation data to row owners. 1555629c3df2SDmitry Karpeev The roots correspond to the dnnz and onnz entries; thus, there are two roots per row. 1556629c3df2SDmitry Karpeev */ 155783b1a929SMark Adams ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 1558629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1559649b366bSFande Kong PetscInt row = bmindices[br], rowowner = 0, brncols, col; 1560629c3df2SDmitry Karpeev const PetscInt *brcols; 1561a4b3d3acSMatthew G Knepley PetscInt rowrel = 0; /* row's relative index on its owner rank */ 1562629c3df2SDmitry Karpeev ierr = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr); 1563649b366bSFande Kong /* how many roots */ 1564649b366bSFande Kong iremote[br].rank = rowowner; iremote[br].index = rowrel; /* edge from bmdnnz to dnnz */ 1565649b366bSFande Kong /* get nonzero pattern */ 156683b1a929SMark Adams ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr); 1567629c3df2SDmitry Karpeev for (k=0; k<brncols; k++) { 1568629c3df2SDmitry Karpeev col = bNindices[brcols[k]]; 1569649b366bSFande Kong if(col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]){ 1570649b366bSFande Kong sub_dnnz[br]++; 1571649b366bSFande Kong }else{ 1572649b366bSFande Kong sub_onnz[br]++; 1573649b366bSFande Kong } 1574629c3df2SDmitry Karpeev } 157583b1a929SMark Adams ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr); 1576629c3df2SDmitry Karpeev } 1577629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1578629c3df2SDmitry Karpeev /* bsf will have to take care of disposing of bedges. */ 1579649b366bSFande Kong ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1580649b366bSFande Kong ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr); 1581649b366bSFande Kong ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr); 1582649b366bSFande Kong ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr); 1583649b366bSFande Kong ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr); 1584649b366bSFande Kong ierr = PetscFree(sub_dnnz);CHKERRQ(ierr); 1585649b366bSFande Kong ierr = PetscFree(sub_onnz);CHKERRQ(ierr); 1586629c3df2SDmitry Karpeev ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr); 1587629c3df2SDmitry Karpeev } 158822d28d08SBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1589629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1590629c3df2SDmitry Karpeev } 1591629c3df2SDmitry Karpeev ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr); 1592629c3df2SDmitry Karpeev ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr); 1593629c3df2SDmitry Karpeev ierr = PetscFree(dnnz);CHKERRQ(ierr); 1594629c3df2SDmitry Karpeev 1595629c3df2SDmitry Karpeev /* Fill by row */ 1596629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1597629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1598629c3df2SDmitry Karpeev IS bNis; 1599629c3df2SDmitry Karpeev PetscInt bN; 1600629c3df2SDmitry Karpeev const PetscInt *bNindices; 1601629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1602629c3df2SDmitry Karpeev ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr); 1603629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1604629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1605629c3df2SDmitry Karpeev Mat B; 1606629c3df2SDmitry Karpeev PetscInt bm, br; 1607629c3df2SDmitry Karpeev const PetscInt *bmindices; 1608629c3df2SDmitry Karpeev B = nest->m[i][j]; 1609629c3df2SDmitry Karpeev if (!B) continue; 1610629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1611629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 161283b1a929SMark Adams ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 1613629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 1614629c3df2SDmitry Karpeev PetscInt row = bmindices[br], brncols, *cols; 1615629c3df2SDmitry Karpeev const PetscInt *brcols; 1616629c3df2SDmitry Karpeev const PetscScalar *brcoldata; 161783b1a929SMark Adams ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1618785e854fSJed Brown ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr); 161926fbe8dcSKarl Rupp for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]]; 1620629c3df2SDmitry Karpeev /* 1621629c3df2SDmitry Karpeev Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match. 1622629c3df2SDmitry Karpeev Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES. 1623629c3df2SDmitry Karpeev */ 1624a2ea699eSBarry Smith ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr); 162583b1a929SMark Adams ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 1626629c3df2SDmitry Karpeev ierr = PetscFree(cols);CHKERRQ(ierr); 1627629c3df2SDmitry Karpeev } 1628629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1629629c3df2SDmitry Karpeev } 1630a2ea699eSBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 1631629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 1632629c3df2SDmitry Karpeev } 1633629c3df2SDmitry Karpeev ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1634629c3df2SDmitry Karpeev ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1635629c3df2SDmitry Karpeev PetscFunctionReturn(0); 1636629c3df2SDmitry Karpeev } 1637629c3df2SDmitry Karpeev 1638659c6bb0SJed Brown /*MC 1639659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 1640659c6bb0SJed Brown 1641659c6bb0SJed Brown Level: intermediate 1642659c6bb0SJed Brown 1643659c6bb0SJed Brown Notes: 1644659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 1645659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 1646950540a4SJed Brown It is usually used with DMComposite and DMCreateMatrix() 1647659c6bb0SJed Brown 1648659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest() 1649659c6bb0SJed Brown M*/ 1650c8883902SJed Brown #undef __FUNCT__ 1651c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest" 16528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A) 1653c8883902SJed Brown { 1654c8883902SJed Brown Mat_Nest *s; 1655c8883902SJed Brown PetscErrorCode ierr; 1656c8883902SJed Brown 1657c8883902SJed Brown PetscFunctionBegin; 1658b00a9115SJed Brown ierr = PetscNewLog(A,&s);CHKERRQ(ierr); 1659c8883902SJed Brown A->data = (void*)s; 1660e7c19651SJed Brown 1661e7c19651SJed Brown s->nr = -1; 1662e7c19651SJed Brown s->nc = -1; 16630298fd71SBarry Smith s->m = NULL; 1664e7c19651SJed Brown s->splitassembly = PETSC_FALSE; 1665c8883902SJed Brown 1666c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 166726fbe8dcSKarl Rupp 1668c8883902SJed Brown A->ops->mult = MatMult_Nest; 16699194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 1670c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 16719194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 1672f8170845SAlex Fikl A->ops->transpose = MatTranspose_Nest; 1673c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 1674c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 1675c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 1676c222c20dSDavid Ham A->ops->copy = MatCopy_Nest; 1677c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 1678c8883902SJed Brown A->ops->getsubmatrix = MatGetSubMatrix_Nest; 1679c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 1680c8883902SJed Brown A->ops->view = MatView_Nest; 1681c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 1682c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 1683c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 1684429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 1685429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 1686a061e289SJed Brown A->ops->scale = MatScale_Nest; 1687a061e289SJed Brown A->ops->shift = MatShift_Nest; 168813135bc6SAlex Fikl A->ops->diagonalset = MatDiagonalSet_Nest; 1689f8170845SAlex Fikl A->ops->setrandom = MatSetRandom_Nest; 1690c8883902SJed Brown 1691c8883902SJed Brown A->spptr = 0; 1692c8883902SJed Brown A->assembled = PETSC_FALSE; 1693c8883902SJed Brown 1694c8883902SJed Brown /* expose Nest api's */ 1695bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C", MatNestGetSubMat_Nest);CHKERRQ(ierr); 1696bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C", MatNestSetSubMat_Nest);CHKERRQ(ierr); 1697bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C", MatNestGetSubMats_Nest);CHKERRQ(ierr); 1698bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C", MatNestGetSize_Nest);CHKERRQ(ierr); 1699bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C", MatNestGetISs_Nest);CHKERRQ(ierr); 1700bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 1701bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C", MatNestSetVecType_Nest);CHKERRQ(ierr); 1702bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C", MatNestSetSubMats_Nest);CHKERRQ(ierr); 170383b1a929SMark Adams ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr); 1704c8883902SJed Brown 1705c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 1706c8883902SJed Brown PetscFunctionReturn(0); 1707c8883902SJed Brown } 1708