1aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I "petscmat.h" I*/ 2b68353e5Sstefano_zampini #include <../src/mat/impls/aij/seq/aij.h> 30c312b8eSJed Brown #include <petscsf.h> 4d8588912SDave May 5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]); 606a1af2fSStefano Zampini static PetscErrorCode MatCreateVecs_Nest(Mat,Vec*,Vec*); 706a1af2fSStefano Zampini static PetscErrorCode MatReset_Nest(Mat); 806a1af2fSStefano Zampini 95e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*); 10c8883902SJed Brown 11d8588912SDave May /* private functions */ 128188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N) 13d8588912SDave May { 14d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 158188e55aSJed Brown PetscInt i,j; 16d8588912SDave May PetscErrorCode ierr; 17d8588912SDave May 18d8588912SDave May PetscFunctionBegin; 198188e55aSJed Brown *m = *n = *M = *N = 0; 208188e55aSJed Brown for (i=0; i<bA->nr; i++) { /* rows */ 218188e55aSJed Brown PetscInt sm,sM; 228188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr); 238188e55aSJed Brown ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr); 248188e55aSJed Brown *m += sm; 258188e55aSJed Brown *M += sM; 26d8588912SDave May } 278188e55aSJed Brown for (j=0; j<bA->nc; j++) { /* cols */ 288188e55aSJed Brown PetscInt sn,sN; 298188e55aSJed Brown ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr); 308188e55aSJed Brown ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr); 318188e55aSJed Brown *n += sn; 328188e55aSJed Brown *N += sN; 33d8588912SDave May } 34d8588912SDave May PetscFunctionReturn(0); 35d8588912SDave May } 36d8588912SDave May 37d8588912SDave May /* operations */ 38207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y) 39d8588912SDave May { 40d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 41207556f9SJed Brown Vec *bx = bA->right,*by = bA->left; 42207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 43d8588912SDave May PetscErrorCode ierr; 44d8588912SDave May 45d8588912SDave May PetscFunctionBegin; 46207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 47207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 48207556f9SJed Brown for (i=0; i<nr; i++) { 49d8588912SDave May ierr = VecZeroEntries(by[i]);CHKERRQ(ierr); 50207556f9SJed Brown for (j=0; j<nc; j++) { 51207556f9SJed Brown if (!bA->m[i][j]) continue; 52d8588912SDave May /* y[i] <- y[i] + A[i][j] * x[j] */ 53d8588912SDave May ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr); 54d8588912SDave May } 55d8588912SDave May } 56207556f9SJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);} 57207556f9SJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 58d8588912SDave May PetscFunctionReturn(0); 59d8588912SDave May } 60d8588912SDave May 619194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z) 629194d70fSJed Brown { 639194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 649194d70fSJed Brown Vec *bx = bA->right,*bz = bA->left; 659194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 669194d70fSJed Brown PetscErrorCode ierr; 679194d70fSJed Brown 689194d70fSJed Brown PetscFunctionBegin; 699194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 709194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 719194d70fSJed Brown for (i=0; i<nr; i++) { 729194d70fSJed Brown if (y != z) { 739194d70fSJed Brown Vec by; 749194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 759194d70fSJed Brown ierr = VecCopy(by,bz[i]);CHKERRQ(ierr); 76336d21e7SJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr); 779194d70fSJed Brown } 789194d70fSJed Brown for (j=0; j<nc; j++) { 799194d70fSJed Brown if (!bA->m[i][j]) continue; 809194d70fSJed Brown /* y[i] <- y[i] + A[i][j] * x[j] */ 819194d70fSJed Brown ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr); 829194d70fSJed Brown } 839194d70fSJed Brown } 849194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);} 859194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);} 869194d70fSJed Brown PetscFunctionReturn(0); 879194d70fSJed Brown } 889194d70fSJed Brown 8952c5f739Sprj- typedef struct { 9052c5f739Sprj- Mat *workC; /* array of Mat with specific containers depending on the underlying MatMatMult implementation */ 9152c5f739Sprj- PetscScalar *tarray; /* buffer for storing all temporary products A[i][j] B[j] */ 9252c5f739Sprj- PetscInt *dm,*dn,k; /* displacements and number of submatrices */ 9352c5f739Sprj- } Nest_Dense; 9452c5f739Sprj- 9552c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat A,Mat B,Mat C) 9652c5f739Sprj- { 9752c5f739Sprj- Mat_Nest *bA = (Mat_Nest*)A->data; 9852c5f739Sprj- PetscContainer container; 9952c5f739Sprj- Nest_Dense *contents; 1004222ddf1SHong Zhang Mat viewB,viewC,seq,productB,workC; 10152c5f739Sprj- const PetscScalar *barray; 10252c5f739Sprj- PetscScalar *carray; 10352c5f739Sprj- PetscInt i,j,M,N,nr = bA->nr,nc = bA->nc,ldb,ldc; 10452c5f739Sprj- PetscErrorCode ierr; 10552c5f739Sprj- 10652c5f739Sprj- PetscFunctionBegin; 10752c5f739Sprj- ierr = PetscObjectQuery((PetscObject)C,"workC",(PetscObject*)&container);CHKERRQ(ierr); 10852c5f739Sprj- if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exist"); 10952c5f739Sprj- ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 11052c5f739Sprj- ierr = MatDenseGetLDA(B,&ldb);CHKERRQ(ierr); 11152c5f739Sprj- ierr = MatDenseGetLDA(C,&ldc);CHKERRQ(ierr); 11252c5f739Sprj- ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr); 11352c5f739Sprj- ierr = MatZeroEntries(C);CHKERRQ(ierr); 11452c5f739Sprj- ierr = MatDenseGetArrayRead(B,&barray);CHKERRQ(ierr); 11552c5f739Sprj- ierr = MatDenseGetArray(C,&carray);CHKERRQ(ierr); 11652c5f739Sprj- for (i=0; i<nr; i++) { 11752c5f739Sprj- ierr = ISGetSize(bA->isglobal.row[i],&M);CHKERRQ(ierr); 11852c5f739Sprj- ierr = MatCreateDense(PetscObjectComm((PetscObject)A),contents->dm[i+1]-contents->dm[i],PETSC_DECIDE,M,N,carray+contents->dm[i],&viewC);CHKERRQ(ierr); 11952c5f739Sprj- ierr = MatDenseGetLocalMatrix(viewC,&seq);CHKERRQ(ierr); 12052c5f739Sprj- ierr = MatSeqDenseSetLDA(seq,ldc);CHKERRQ(ierr); 12152c5f739Sprj- for (j=0; j<nc; j++) { 12252c5f739Sprj- if (!bA->m[i][j]) continue; 12352c5f739Sprj- ierr = ISGetSize(bA->isglobal.col[j],&M);CHKERRQ(ierr); 12452c5f739Sprj- ierr = MatCreateDense(PetscObjectComm((PetscObject)A),contents->dn[j+1]-contents->dn[j],PETSC_DECIDE,M,N,(PetscScalar*)(barray+contents->dn[j]),&viewB);CHKERRQ(ierr); 12552c5f739Sprj- ierr = MatDenseGetLocalMatrix(viewB,&seq);CHKERRQ(ierr); 12652c5f739Sprj- ierr = MatSeqDenseSetLDA(seq,ldb);CHKERRQ(ierr); 1274222ddf1SHong Zhang 1284222ddf1SHong Zhang /* MatMatMultNumeric(bA->m[i][j],viewB,contents->workC[i*nc + j]); */ 1294222ddf1SHong Zhang workC = contents->workC[i*nc + j]; 1304222ddf1SHong Zhang productB = workC->product->B; 1314222ddf1SHong Zhang workC->product->B = viewB; /* use newly created dense matrix viewB */ 132*18992e5dSStefano Zampini ierr = (*workC->ops->productnumeric)(workC);CHKERRQ(ierr); 13352c5f739Sprj- ierr = MatDestroy(&viewB);CHKERRQ(ierr); 1344222ddf1SHong Zhang workC->product->B = productB; /* resume original B */ 1354222ddf1SHong Zhang 13652c5f739Sprj- /* C[i] <- workC + C[i] */ 13752c5f739Sprj- ierr = MatAXPY(viewC,1.0,contents->workC[i*nc + j],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 13852c5f739Sprj- } 13952c5f739Sprj- ierr = MatDestroy(&viewC);CHKERRQ(ierr); 14052c5f739Sprj- } 14152c5f739Sprj- ierr = MatDenseRestoreArray(C,&carray);CHKERRQ(ierr); 14252c5f739Sprj- ierr = MatDenseRestoreArrayRead(B,&barray);CHKERRQ(ierr); 1434222ddf1SHong Zhang 1444222ddf1SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1454222ddf1SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14652c5f739Sprj- PetscFunctionReturn(0); 14752c5f739Sprj- } 14852c5f739Sprj- 14952c5f739Sprj- PetscErrorCode MatNest_DenseDestroy(void *ctx) 15052c5f739Sprj- { 15152c5f739Sprj- Nest_Dense *contents = (Nest_Dense*)ctx; 15252c5f739Sprj- PetscInt i; 15352c5f739Sprj- PetscErrorCode ierr; 15452c5f739Sprj- 15552c5f739Sprj- PetscFunctionBegin; 15652c5f739Sprj- ierr = PetscFree(contents->tarray);CHKERRQ(ierr); 15752c5f739Sprj- for (i=0; i<contents->k; i++) { 15852c5f739Sprj- ierr = MatDestroy(contents->workC + i);CHKERRQ(ierr); 15952c5f739Sprj- } 16052c5f739Sprj- ierr = PetscFree3(contents->dm,contents->dn,contents->workC);CHKERRQ(ierr); 16152c5f739Sprj- ierr = PetscFree(contents);CHKERRQ(ierr); 16252c5f739Sprj- PetscFunctionReturn(0); 16352c5f739Sprj- } 16452c5f739Sprj- 1654222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat A,Mat B,PetscReal fill,Mat C) 16652c5f739Sprj- { 16752c5f739Sprj- Mat_Nest *bA = (Mat_Nest*)A->data; 1684222ddf1SHong Zhang Mat viewB,viewSeq,workC; 16952c5f739Sprj- const PetscScalar *barray; 17052c5f739Sprj- PetscInt i,j,M,N,m,nr = bA->nr,nc = bA->nc,maxm = 0,ldb; 17152c5f739Sprj- PetscContainer container; 1724222ddf1SHong Zhang Nest_Dense *contents=NULL; 17352c5f739Sprj- PetscErrorCode ierr; 17452c5f739Sprj- 17552c5f739Sprj- PetscFunctionBegin; 17652c5f739Sprj- ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr); 1771a4d94e9SStefano Zampini if (!C->assembled) { 17852c5f739Sprj- ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 17952c5f739Sprj- ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr); 1804222ddf1SHong Zhang 1814222ddf1SHong Zhang ierr = MatSetSizes(C,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1824222ddf1SHong Zhang ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 183*18992e5dSStefano Zampini ierr = MatSetUp(C);CHKERRQ(ierr); 18452c5f739Sprj- } 18552c5f739Sprj- 18652c5f739Sprj- ierr = PetscNew(&contents);CHKERRQ(ierr); 18752c5f739Sprj- ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 18852c5f739Sprj- ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 18952c5f739Sprj- ierr = PetscContainerSetUserDestroy(container,MatNest_DenseDestroy);CHKERRQ(ierr); 1904222ddf1SHong Zhang ierr = PetscObjectCompose((PetscObject)C,"workC",(PetscObject)container);CHKERRQ(ierr); 19152c5f739Sprj- ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 19252c5f739Sprj- ierr = PetscCalloc3(nr+1,&contents->dm,nc+1,&contents->dn,nr*nc,&contents->workC);CHKERRQ(ierr); 19352c5f739Sprj- contents->k = nr*nc; 19452c5f739Sprj- for (i=0; i<nr; i++) { 19552c5f739Sprj- ierr = ISGetLocalSize(bA->isglobal.row[i],contents->dm + i+1);CHKERRQ(ierr); 19652c5f739Sprj- maxm = PetscMax(maxm,contents->dm[i+1]); 19752c5f739Sprj- contents->dm[i+1] += contents->dm[i]; 19852c5f739Sprj- } 19952c5f739Sprj- for (i=0; i<nc; i++) { 20052c5f739Sprj- ierr = ISGetLocalSize(bA->isglobal.col[i],contents->dn + i+1);CHKERRQ(ierr); 20152c5f739Sprj- contents->dn[i+1] += contents->dn[i]; 20252c5f739Sprj- } 20352c5f739Sprj- ierr = PetscMalloc1(maxm*N,&contents->tarray);CHKERRQ(ierr); 20452c5f739Sprj- ierr = MatDenseGetLDA(B,&ldb);CHKERRQ(ierr); 20552c5f739Sprj- ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr); 20652c5f739Sprj- ierr = MatDenseGetArrayRead(B,&barray);CHKERRQ(ierr); 20752c5f739Sprj- /* loops are permuted compared to MatMatMultNumeric so that viewB is created only once per column of A */ 20852c5f739Sprj- for (j=0; j<nc; j++) { 20952c5f739Sprj- ierr = ISGetSize(bA->isglobal.col[j],&M);CHKERRQ(ierr); 21052c5f739Sprj- ierr = MatCreateDense(PetscObjectComm((PetscObject)A),contents->dn[j+1]-contents->dn[j],PETSC_DECIDE,M,N,(PetscScalar*)(barray+contents->dn[j]),&viewB);CHKERRQ(ierr); 21152c5f739Sprj- ierr = MatDenseGetLocalMatrix(viewB,&viewSeq);CHKERRQ(ierr); 21252c5f739Sprj- ierr = MatSeqDenseSetLDA(viewSeq,ldb);CHKERRQ(ierr); 21352c5f739Sprj- for (i=0; i<nr; i++) { 21452c5f739Sprj- if (!bA->m[i][j]) continue; 21552c5f739Sprj- /* MatMatMultSymbolic may attach a specific container (depending on MatType of bA->m[i][j]) to workC[i][j] */ 2164222ddf1SHong Zhang 2174222ddf1SHong Zhang ierr = MatProductCreate(bA->m[i][j],viewB,NULL,&contents->workC[i*nc + j]);CHKERRQ(ierr); 2184222ddf1SHong Zhang workC = contents->workC[i*nc + j]; 2194222ddf1SHong Zhang ierr = MatProductSetType(workC,MATPRODUCT_AB);CHKERRQ(ierr); 2204222ddf1SHong Zhang ierr = MatProductSetAlgorithm(workC,"default");CHKERRQ(ierr); 2214222ddf1SHong Zhang ierr = MatProductSetFill(workC,fill);CHKERRQ(ierr); 2224222ddf1SHong Zhang ierr = MatProductSetFromOptions(workC);CHKERRQ(ierr); 2234222ddf1SHong Zhang ierr = MatProductSymbolic(workC);CHKERRQ(ierr); 2244222ddf1SHong Zhang 2254222ddf1SHong Zhang ierr = MatDenseGetLocalMatrix(workC,&viewSeq);CHKERRQ(ierr); 22652c5f739Sprj- /* free the memory allocated in MatMatMultSymbolic, since tarray will be shared by all Mat */ 22752c5f739Sprj- ierr = MatSeqDenseSetPreallocation(viewSeq,contents->tarray);CHKERRQ(ierr); 22852c5f739Sprj- } 22952c5f739Sprj- ierr = MatDestroy(&viewB);CHKERRQ(ierr); 23052c5f739Sprj- } 23152c5f739Sprj- ierr = MatDenseRestoreArrayRead(B,&barray);CHKERRQ(ierr); 23252c5f739Sprj- 2334222ddf1SHong Zhang C->ops->matmultnumeric = MatMatMultNumeric_Nest_Dense; 23452c5f739Sprj- PetscFunctionReturn(0); 23552c5f739Sprj- } 23652c5f739Sprj- 2374222ddf1SHong Zhang /* --------------------------------------------------------- */ 2384222ddf1SHong Zhang static PetscErrorCode MatProductSetFromOptions_Nest_Dense_AB(Mat C) 2394222ddf1SHong Zhang { 2404222ddf1SHong Zhang PetscFunctionBegin; 2414222ddf1SHong Zhang C->ops->matmultsymbolic = MatMatMultSymbolic_Nest_Dense; 2424222ddf1SHong Zhang C->ops->productsymbolic = MatProductSymbolic_AB; 2434222ddf1SHong Zhang PetscFunctionReturn(0); 2444222ddf1SHong Zhang } 2454222ddf1SHong Zhang 2464222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Nest_Dense(Mat C) 24752c5f739Sprj- { 24852c5f739Sprj- PetscErrorCode ierr; 2494222ddf1SHong Zhang Mat_Product *product = C->product; 25052c5f739Sprj- 25152c5f739Sprj- PetscFunctionBegin; 2524222ddf1SHong Zhang if (product->type == MATPRODUCT_AB) { 2534222ddf1SHong Zhang ierr = MatProductSetFromOptions_Nest_Dense_AB(C);CHKERRQ(ierr); 254544a5e07SHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"MatProduct type %s is not supported for Nest and Dense matrices",MatProductTypes[product->type]); 25552c5f739Sprj- PetscFunctionReturn(0); 25652c5f739Sprj- } 2574222ddf1SHong Zhang /* --------------------------------------------------------- */ 25852c5f739Sprj- 259207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y) 260d8588912SDave May { 261d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 262207556f9SJed Brown Vec *bx = bA->left,*by = bA->right; 263207556f9SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 264d8588912SDave May PetscErrorCode ierr; 265d8588912SDave May 266d8588912SDave May PetscFunctionBegin; 267609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 268609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 269207556f9SJed Brown for (j=0; j<nc; j++) { 270609e31cbSJed Brown ierr = VecZeroEntries(by[j]);CHKERRQ(ierr); 271609e31cbSJed Brown for (i=0; i<nr; i++) { 2726c75ac25SJed Brown if (!bA->m[i][j]) continue; 273609e31cbSJed Brown /* y[j] <- y[j] + (A[i][j])^T * x[i] */ 274609e31cbSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr); 275d8588912SDave May } 276d8588912SDave May } 277609e31cbSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 278609e31cbSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);} 279d8588912SDave May PetscFunctionReturn(0); 280d8588912SDave May } 281d8588912SDave May 2829194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z) 2839194d70fSJed Brown { 2849194d70fSJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 2859194d70fSJed Brown Vec *bx = bA->left,*bz = bA->right; 2869194d70fSJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 2879194d70fSJed Brown PetscErrorCode ierr; 2889194d70fSJed Brown 2899194d70fSJed Brown PetscFunctionBegin; 2909194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 2919194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 2929194d70fSJed Brown for (j=0; j<nc; j++) { 2939194d70fSJed Brown if (y != z) { 2949194d70fSJed Brown Vec by; 2959194d70fSJed Brown ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 2969194d70fSJed Brown ierr = VecCopy(by,bz[j]);CHKERRQ(ierr); 2979194d70fSJed Brown ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr); 2989194d70fSJed Brown } 2999194d70fSJed Brown for (i=0; i<nr; i++) { 3006c75ac25SJed Brown if (!bA->m[i][j]) continue; 3019194d70fSJed Brown /* z[j] <- y[j] + (A[i][j])^T * x[i] */ 3029194d70fSJed Brown ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr); 3039194d70fSJed Brown } 3049194d70fSJed Brown } 3059194d70fSJed Brown for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);} 3069194d70fSJed Brown for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);} 3079194d70fSJed Brown PetscFunctionReturn(0); 3089194d70fSJed Brown } 3099194d70fSJed Brown 310f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B) 311f8170845SAlex Fikl { 312f8170845SAlex Fikl Mat_Nest *bA = (Mat_Nest*)A->data, *bC; 313f8170845SAlex Fikl Mat C; 314f8170845SAlex Fikl PetscInt i,j,nr = bA->nr,nc = bA->nc; 315f8170845SAlex Fikl PetscErrorCode ierr; 316f8170845SAlex Fikl 317f8170845SAlex Fikl PetscFunctionBegin; 318cf37664fSBarry Smith if (reuse == MAT_INPLACE_MATRIX && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place"); 319f8170845SAlex Fikl 320cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 321f8170845SAlex Fikl Mat *subs; 322f8170845SAlex Fikl IS *is_row,*is_col; 323f8170845SAlex Fikl 324f8170845SAlex Fikl ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr); 325f8170845SAlex Fikl ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr); 326f8170845SAlex Fikl ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr); 327cf37664fSBarry Smith if (reuse == MAT_INPLACE_MATRIX) { 328ddeb9bd8SAlex Fikl for (i=0; i<nr; i++) { 329ddeb9bd8SAlex Fikl for (j=0; j<nc; j++) { 330ddeb9bd8SAlex Fikl subs[i + nr * j] = bA->m[i][j]; 331ddeb9bd8SAlex Fikl } 332ddeb9bd8SAlex Fikl } 333ddeb9bd8SAlex Fikl } 334ddeb9bd8SAlex Fikl 335f8170845SAlex Fikl ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr); 336f8170845SAlex Fikl ierr = PetscFree(subs);CHKERRQ(ierr); 3373d994f23SBarry Smith ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr); 338f8170845SAlex Fikl } else { 339f8170845SAlex Fikl C = *B; 340f8170845SAlex Fikl } 341f8170845SAlex Fikl 342f8170845SAlex Fikl bC = (Mat_Nest*)C->data; 343f8170845SAlex Fikl for (i=0; i<nr; i++) { 344f8170845SAlex Fikl for (j=0; j<nc; j++) { 345f8170845SAlex Fikl if (bA->m[i][j]) { 346f8170845SAlex Fikl ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr); 347f8170845SAlex Fikl } else { 348f8170845SAlex Fikl bC->m[j][i] = NULL; 349f8170845SAlex Fikl } 350f8170845SAlex Fikl } 351f8170845SAlex Fikl } 352f8170845SAlex Fikl 353cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 354f8170845SAlex Fikl *B = C; 355f8170845SAlex Fikl } else { 356f8170845SAlex Fikl ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr); 357f8170845SAlex Fikl } 358f8170845SAlex Fikl PetscFunctionReturn(0); 359f8170845SAlex Fikl } 360f8170845SAlex Fikl 361e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list) 362e2d7f03fSJed Brown { 363e2d7f03fSJed Brown PetscErrorCode ierr; 364e2d7f03fSJed Brown IS *lst = *list; 365e2d7f03fSJed Brown PetscInt i; 366e2d7f03fSJed Brown 367e2d7f03fSJed Brown PetscFunctionBegin; 368e2d7f03fSJed Brown if (!lst) PetscFunctionReturn(0); 3696bf464f9SBarry Smith for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);} 370e2d7f03fSJed Brown ierr = PetscFree(lst);CHKERRQ(ierr); 3710298fd71SBarry Smith *list = NULL; 372e2d7f03fSJed Brown PetscFunctionReturn(0); 373e2d7f03fSJed Brown } 374e2d7f03fSJed Brown 37506a1af2fSStefano Zampini static PetscErrorCode MatReset_Nest(Mat A) 376d8588912SDave May { 377d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 378d8588912SDave May PetscInt i,j; 379d8588912SDave May PetscErrorCode ierr; 380d8588912SDave May 381d8588912SDave May PetscFunctionBegin; 382d8588912SDave May /* release the matrices and the place holders */ 383e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr); 384e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr); 385e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 386e2d7f03fSJed Brown ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 387d8588912SDave May 388d8588912SDave May ierr = PetscFree(vs->row_len);CHKERRQ(ierr); 389d8588912SDave May ierr = PetscFree(vs->col_len);CHKERRQ(ierr); 39006a1af2fSStefano Zampini ierr = PetscFree(vs->nnzstate);CHKERRQ(ierr); 391d8588912SDave May 392207556f9SJed Brown ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr); 393207556f9SJed Brown 394d8588912SDave May /* release the matrices and the place holders */ 395d8588912SDave May if (vs->m) { 396d8588912SDave May for (i=0; i<vs->nr; i++) { 397d8588912SDave May for (j=0; j<vs->nc; j++) { 3986bf464f9SBarry Smith ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr); 399d8588912SDave May } 400d8588912SDave May ierr = PetscFree(vs->m[i]);CHKERRQ(ierr); 401d8588912SDave May } 402d8588912SDave May ierr = PetscFree(vs->m);CHKERRQ(ierr); 403d8588912SDave May } 40406a1af2fSStefano Zampini 40506a1af2fSStefano Zampini /* restore defaults */ 40606a1af2fSStefano Zampini vs->nr = 0; 40706a1af2fSStefano Zampini vs->nc = 0; 40806a1af2fSStefano Zampini vs->splitassembly = PETSC_FALSE; 40906a1af2fSStefano Zampini PetscFunctionReturn(0); 41006a1af2fSStefano Zampini } 41106a1af2fSStefano Zampini 41206a1af2fSStefano Zampini static PetscErrorCode MatDestroy_Nest(Mat A) 41306a1af2fSStefano Zampini { 41406a1af2fSStefano Zampini PetscErrorCode ierr; 41506a1af2fSStefano Zampini 41606a1af2fSStefano Zampini ierr = MatReset_Nest(A);CHKERRQ(ierr); 417bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 418d8588912SDave May 419bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr); 420bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr); 421bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr); 422bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr); 423bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr); 424bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr); 425bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr); 426bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr); 4270899c546SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",0);CHKERRQ(ierr); 4280899c546SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",0);CHKERRQ(ierr); 4295e3038f0Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr); 4305e3038f0Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr); 4314222ddf1SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_seqdense_C",NULL);CHKERRQ(ierr); 4324222ddf1SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_mpidense_C",NULL);CHKERRQ(ierr); 4334222ddf1SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_dense_C",NULL);CHKERRQ(ierr); 434d8588912SDave May PetscFunctionReturn(0); 435d8588912SDave May } 436d8588912SDave May 437381b8e50SStefano Zampini static PetscErrorCode MatMissingDiagonal_Nest(Mat mat,PetscBool *missing,PetscInt *dd) 438381b8e50SStefano Zampini { 439381b8e50SStefano Zampini Mat_Nest *vs = (Mat_Nest*)mat->data; 440381b8e50SStefano Zampini PetscInt i; 441381b8e50SStefano Zampini PetscErrorCode ierr; 442381b8e50SStefano Zampini 443381b8e50SStefano Zampini PetscFunctionBegin; 444381b8e50SStefano Zampini if (dd) *dd = 0; 445381b8e50SStefano Zampini if (!vs->nr) { 446381b8e50SStefano Zampini *missing = PETSC_TRUE; 447381b8e50SStefano Zampini PetscFunctionReturn(0); 448381b8e50SStefano Zampini } 449381b8e50SStefano Zampini *missing = PETSC_FALSE; 450381b8e50SStefano Zampini for (i = 0; i < vs->nr && !(*missing); i++) { 451381b8e50SStefano Zampini *missing = PETSC_TRUE; 452381b8e50SStefano Zampini if (vs->m[i][i]) { 453381b8e50SStefano Zampini ierr = MatMissingDiagonal(vs->m[i][i],missing,NULL);CHKERRQ(ierr); 454381b8e50SStefano Zampini if (*missing && dd) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"First missing entry not yet implemented"); 455381b8e50SStefano Zampini } 456381b8e50SStefano Zampini } 457381b8e50SStefano Zampini PetscFunctionReturn(0); 458381b8e50SStefano Zampini } 459381b8e50SStefano Zampini 460207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type) 461d8588912SDave May { 462d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 463d8588912SDave May PetscInt i,j; 464d8588912SDave May PetscErrorCode ierr; 46506a1af2fSStefano Zampini PetscBool nnzstate = PETSC_FALSE; 466d8588912SDave May 467d8588912SDave May PetscFunctionBegin; 468d8588912SDave May for (i=0; i<vs->nr; i++) { 469d8588912SDave May for (j=0; j<vs->nc; j++) { 47006a1af2fSStefano Zampini PetscObjectState subnnzstate = 0; 471e7c19651SJed Brown if (vs->m[i][j]) { 472e7c19651SJed Brown ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); 473e7c19651SJed Brown if (!vs->splitassembly) { 474e7c19651SJed Brown /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested 475e7c19651SJed 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 476e7c19651SJed Brown * already performing an assembly, but the result would by more complicated and appears to offer less 477e7c19651SJed Brown * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an 478e7c19651SJed Brown * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives". 479e7c19651SJed Brown */ 480e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 48106a1af2fSStefano Zampini ierr = MatGetNonzeroState(vs->m[i][j],&subnnzstate);CHKERRQ(ierr); 482e7c19651SJed Brown } 483e7c19651SJed Brown } 48406a1af2fSStefano Zampini nnzstate = (PetscBool)(nnzstate || vs->nnzstate[i*vs->nc+j] != subnnzstate); 48506a1af2fSStefano Zampini vs->nnzstate[i*vs->nc+j] = subnnzstate; 486d8588912SDave May } 487d8588912SDave May } 48806a1af2fSStefano Zampini if (nnzstate) A->nonzerostate++; 489d8588912SDave May PetscFunctionReturn(0); 490d8588912SDave May } 491d8588912SDave May 492207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type) 493d8588912SDave May { 494d8588912SDave May Mat_Nest *vs = (Mat_Nest*)A->data; 495d8588912SDave May PetscInt i,j; 496d8588912SDave May PetscErrorCode ierr; 497d8588912SDave May 498d8588912SDave May PetscFunctionBegin; 499d8588912SDave May for (i=0; i<vs->nr; i++) { 500d8588912SDave May for (j=0; j<vs->nc; j++) { 501e7c19651SJed Brown if (vs->m[i][j]) { 502e7c19651SJed Brown if (vs->splitassembly) { 503e7c19651SJed Brown ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); 504e7c19651SJed Brown } 505e7c19651SJed Brown } 506d8588912SDave May } 507d8588912SDave May } 508d8588912SDave May PetscFunctionReturn(0); 509d8588912SDave May } 510d8588912SDave May 511f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B) 512d8588912SDave May { 513207556f9SJed Brown PetscErrorCode ierr; 514f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 515f349c1fdSJed Brown PetscInt j; 516f349c1fdSJed Brown Mat sub; 517d8588912SDave May 518d8588912SDave May PetscFunctionBegin; 5190298fd71SBarry Smith sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */ 520f349c1fdSJed Brown for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j]; 5214994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 522f349c1fdSJed Brown *B = sub; 523f349c1fdSJed Brown PetscFunctionReturn(0); 524d8588912SDave May } 525d8588912SDave May 526f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B) 527f349c1fdSJed Brown { 528207556f9SJed Brown PetscErrorCode ierr; 529f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 530f349c1fdSJed Brown PetscInt i; 531f349c1fdSJed Brown Mat sub; 532f349c1fdSJed Brown 533f349c1fdSJed Brown PetscFunctionBegin; 5340298fd71SBarry Smith sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */ 535f349c1fdSJed Brown for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col]; 5364994cf47SJed Brown if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */ 537f349c1fdSJed Brown *B = sub; 538f349c1fdSJed Brown PetscFunctionReturn(0); 539d8588912SDave May } 540d8588912SDave May 541f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found) 542f349c1fdSJed Brown { 543f349c1fdSJed Brown PetscErrorCode ierr; 544f349c1fdSJed Brown PetscInt i; 545f349c1fdSJed Brown PetscBool flg; 546f349c1fdSJed Brown 547f349c1fdSJed Brown PetscFunctionBegin; 548f349c1fdSJed Brown PetscValidPointer(list,3); 549f349c1fdSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,4); 550f349c1fdSJed Brown PetscValidIntPointer(found,5); 551f349c1fdSJed Brown *found = -1; 552f349c1fdSJed Brown for (i=0; i<n; i++) { 553207556f9SJed Brown if (!list[i]) continue; 554320466b0SStefano Zampini ierr = ISEqualUnsorted(list[i],is,&flg);CHKERRQ(ierr); 555f349c1fdSJed Brown if (flg) { 556f349c1fdSJed Brown *found = i; 557f349c1fdSJed Brown PetscFunctionReturn(0); 558f349c1fdSJed Brown } 559f349c1fdSJed Brown } 560ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set"); 561f349c1fdSJed Brown PetscFunctionReturn(0); 562f349c1fdSJed Brown } 563f349c1fdSJed Brown 5648188e55aSJed Brown /* Get a block row as a new MatNest */ 5658188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B) 5668188e55aSJed Brown { 5678188e55aSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 5688188e55aSJed Brown char keyname[256]; 5698188e55aSJed Brown PetscErrorCode ierr; 5708188e55aSJed Brown 5718188e55aSJed Brown PetscFunctionBegin; 5720298fd71SBarry Smith *B = NULL; 5738caf3d72SBarry Smith ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr); 5748188e55aSJed Brown ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr); 5758188e55aSJed Brown if (*B) PetscFunctionReturn(0); 5768188e55aSJed Brown 577ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr); 57826fbe8dcSKarl Rupp 5798188e55aSJed Brown (*B)->assembled = A->assembled; 58026fbe8dcSKarl Rupp 5818188e55aSJed Brown ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr); 5828188e55aSJed Brown ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */ 5838188e55aSJed Brown PetscFunctionReturn(0); 5848188e55aSJed Brown } 5858188e55aSJed Brown 586f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B) 587f349c1fdSJed Brown { 588f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 5898188e55aSJed Brown PetscErrorCode ierr; 5906b3a5b13SJed Brown PetscInt row,col; 591e072481dSJed Brown PetscBool same,isFullCol,isFullColGlobal; 592f349c1fdSJed Brown 593f349c1fdSJed Brown PetscFunctionBegin; 5948188e55aSJed Brown /* Check if full column space. This is a hack */ 5958188e55aSJed Brown isFullCol = PETSC_FALSE; 596251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr); 5978188e55aSJed Brown if (same) { 59877019fcaSJed Brown PetscInt n,first,step,i,an,am,afirst,astep; 5998188e55aSJed Brown ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 6008188e55aSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 60177019fcaSJed Brown isFullCol = PETSC_TRUE; 60205ce4453SJed Brown for (i=0,an=A->cmap->rstart; i<vs->nc; i++) { 60306a1af2fSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)is->col[i],ISSTRIDE,&same);CHKERRQ(ierr); 60477019fcaSJed Brown ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr); 60506a1af2fSStefano Zampini if (same) { 60606a1af2fSStefano Zampini ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr); 60777019fcaSJed Brown if (afirst != an || astep != step) isFullCol = PETSC_FALSE; 60806a1af2fSStefano Zampini } else isFullCol = PETSC_FALSE; 60977019fcaSJed Brown an += am; 61077019fcaSJed Brown } 61105ce4453SJed Brown if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE; 6128188e55aSJed Brown } 613b2566f29SBarry Smith ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr); 6148188e55aSJed Brown 615427230ceSLisandro Dalcin if (isFullColGlobal && vs->nc > 1) { 6168188e55aSJed Brown PetscInt row; 6178188e55aSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 6188188e55aSJed Brown ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr); 6198188e55aSJed Brown } else { 620f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr); 621f349c1fdSJed Brown ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr); 622b6480e04SStefano Zampini if (!vs->m[row][col]) { 623b6480e04SStefano Zampini PetscInt lr,lc; 624b6480e04SStefano Zampini 625b6480e04SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr); 626b6480e04SStefano Zampini ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr); 627b6480e04SStefano Zampini ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr); 628b6480e04SStefano Zampini ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 629fa9f0909SStefano Zampini ierr = MatSetType(vs->m[row][col],MATAIJ);CHKERRQ(ierr); 630fa9f0909SStefano Zampini ierr = MatSeqAIJSetPreallocation(vs->m[row][col],0,NULL);CHKERRQ(ierr); 631fa9f0909SStefano Zampini ierr = MatMPIAIJSetPreallocation(vs->m[row][col],0,NULL,0,NULL);CHKERRQ(ierr); 632b6480e04SStefano Zampini ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr); 633b6480e04SStefano Zampini ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 634b6480e04SStefano Zampini ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 635b6480e04SStefano Zampini } 636f349c1fdSJed Brown *B = vs->m[row][col]; 6378188e55aSJed Brown } 638f349c1fdSJed Brown PetscFunctionReturn(0); 639f349c1fdSJed Brown } 640f349c1fdSJed Brown 64106a1af2fSStefano Zampini /* 64206a1af2fSStefano Zampini TODO: This does not actually returns a submatrix we can modify 64306a1af2fSStefano Zampini */ 6447dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B) 645f349c1fdSJed Brown { 646f349c1fdSJed Brown PetscErrorCode ierr; 647f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 648f349c1fdSJed Brown Mat sub; 649f349c1fdSJed Brown 650f349c1fdSJed Brown PetscFunctionBegin; 651f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr); 652f349c1fdSJed Brown switch (reuse) { 653f349c1fdSJed Brown case MAT_INITIAL_MATRIX: 6547874fa86SDave May if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); } 655f349c1fdSJed Brown *B = sub; 656f349c1fdSJed Brown break; 657f349c1fdSJed Brown case MAT_REUSE_MATRIX: 658ce94432eSBarry Smith if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call"); 659f349c1fdSJed Brown break; 660f349c1fdSJed Brown case MAT_IGNORE_MATRIX: /* Nothing to do */ 661f349c1fdSJed Brown break; 662511c6705SHong Zhang case MAT_INPLACE_MATRIX: /* Nothing to do */ 663511c6705SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet"); 664511c6705SHong Zhang break; 665f349c1fdSJed Brown } 666f349c1fdSJed Brown PetscFunctionReturn(0); 667f349c1fdSJed Brown } 668f349c1fdSJed Brown 669f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 670f349c1fdSJed Brown { 671f349c1fdSJed Brown PetscErrorCode ierr; 672f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 673f349c1fdSJed Brown Mat sub; 674f349c1fdSJed Brown 675f349c1fdSJed Brown PetscFunctionBegin; 676f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 677f349c1fdSJed Brown /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */ 678f349c1fdSJed Brown if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);} 679f349c1fdSJed Brown *B = sub; 680d8588912SDave May PetscFunctionReturn(0); 681d8588912SDave May } 682d8588912SDave May 683207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B) 684d8588912SDave May { 685d8588912SDave May PetscErrorCode ierr; 686f349c1fdSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 687f349c1fdSJed Brown Mat sub; 688d8588912SDave May 689d8588912SDave May PetscFunctionBegin; 690f349c1fdSJed Brown ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr); 691ce94432eSBarry Smith if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten"); 692f349c1fdSJed Brown if (sub) { 693ce94432eSBarry Smith if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times"); 6946bf464f9SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 695d8588912SDave May } 696d8588912SDave May PetscFunctionReturn(0); 697d8588912SDave May } 698d8588912SDave May 6997874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v) 7007874fa86SDave May { 7017874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 7027874fa86SDave May PetscInt i; 7037874fa86SDave May PetscErrorCode ierr; 7047874fa86SDave May 7057874fa86SDave May PetscFunctionBegin; 7067874fa86SDave May for (i=0; i<bA->nr; i++) { 707429bac76SJed Brown Vec bv; 708429bac76SJed Brown ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 7097874fa86SDave May if (bA->m[i][i]) { 710429bac76SJed Brown ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr); 7117874fa86SDave May } else { 7125159a857SMatthew G. Knepley ierr = VecSet(bv,0.0);CHKERRQ(ierr); 7137874fa86SDave May } 714429bac76SJed Brown ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 7157874fa86SDave May } 7167874fa86SDave May PetscFunctionReturn(0); 7177874fa86SDave May } 7187874fa86SDave May 7197874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r) 7207874fa86SDave May { 7217874fa86SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 722429bac76SJed Brown Vec bl,*br; 7237874fa86SDave May PetscInt i,j; 7247874fa86SDave May PetscErrorCode ierr; 7257874fa86SDave May 7267874fa86SDave May PetscFunctionBegin; 7273f800ebeSJed Brown ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr); 7282e6472ebSElliott Sales de Andrade if (r) { 729429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 7302e6472ebSElliott Sales de Andrade } 7312e6472ebSElliott Sales de Andrade bl = NULL; 7327874fa86SDave May for (i=0; i<bA->nr; i++) { 7332e6472ebSElliott Sales de Andrade if (l) { 734429bac76SJed Brown ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 7352e6472ebSElliott Sales de Andrade } 7367874fa86SDave May for (j=0; j<bA->nc; j++) { 7377874fa86SDave May if (bA->m[i][j]) { 738429bac76SJed Brown ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr); 7397874fa86SDave May } 7407874fa86SDave May } 7412e6472ebSElliott Sales de Andrade if (l) { 742a061e289SJed Brown ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr); 7437874fa86SDave May } 7442e6472ebSElliott Sales de Andrade } 7452e6472ebSElliott Sales de Andrade if (r) { 746429bac76SJed Brown for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);} 7472e6472ebSElliott Sales de Andrade } 748429bac76SJed Brown ierr = PetscFree(br);CHKERRQ(ierr); 7497874fa86SDave May PetscFunctionReturn(0); 7507874fa86SDave May } 7517874fa86SDave May 752a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a) 753a061e289SJed Brown { 754a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 755a061e289SJed Brown PetscInt i,j; 756a061e289SJed Brown PetscErrorCode ierr; 757a061e289SJed Brown 758a061e289SJed Brown PetscFunctionBegin; 759a061e289SJed Brown for (i=0; i<bA->nr; i++) { 760a061e289SJed Brown for (j=0; j<bA->nc; j++) { 761a061e289SJed Brown if (bA->m[i][j]) { 762a061e289SJed Brown ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr); 763a061e289SJed Brown } 764a061e289SJed Brown } 765a061e289SJed Brown } 766a061e289SJed Brown PetscFunctionReturn(0); 767a061e289SJed Brown } 768a061e289SJed Brown 769a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a) 770a061e289SJed Brown { 771a061e289SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 772a061e289SJed Brown PetscInt i; 773a061e289SJed Brown PetscErrorCode ierr; 77406a1af2fSStefano Zampini PetscBool nnzstate = PETSC_FALSE; 775a061e289SJed Brown 776a061e289SJed Brown PetscFunctionBegin; 777a061e289SJed Brown for (i=0; i<bA->nr; i++) { 77806a1af2fSStefano Zampini PetscObjectState subnnzstate = 0; 779ce94432eSBarry 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); 780a061e289SJed Brown ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr); 78106a1af2fSStefano Zampini ierr = MatGetNonzeroState(bA->m[i][i],&subnnzstate);CHKERRQ(ierr); 78206a1af2fSStefano Zampini nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate); 78306a1af2fSStefano Zampini bA->nnzstate[i*bA->nc+i] = subnnzstate; 784a061e289SJed Brown } 78506a1af2fSStefano Zampini if (nnzstate) A->nonzerostate++; 786a061e289SJed Brown PetscFunctionReturn(0); 787a061e289SJed Brown } 788a061e289SJed Brown 78913135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is) 79013135bc6SAlex Fikl { 79113135bc6SAlex Fikl Mat_Nest *bA = (Mat_Nest*)A->data; 79213135bc6SAlex Fikl PetscInt i; 79313135bc6SAlex Fikl PetscErrorCode ierr; 79406a1af2fSStefano Zampini PetscBool nnzstate = PETSC_FALSE; 79513135bc6SAlex Fikl 79613135bc6SAlex Fikl PetscFunctionBegin; 79713135bc6SAlex Fikl for (i=0; i<bA->nr; i++) { 79806a1af2fSStefano Zampini PetscObjectState subnnzstate = 0; 79913135bc6SAlex Fikl Vec bv; 80013135bc6SAlex Fikl ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 80113135bc6SAlex Fikl if (bA->m[i][i]) { 80213135bc6SAlex Fikl ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr); 80306a1af2fSStefano Zampini ierr = MatGetNonzeroState(bA->m[i][i],&subnnzstate);CHKERRQ(ierr); 80413135bc6SAlex Fikl } 80513135bc6SAlex Fikl ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr); 80606a1af2fSStefano Zampini nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate); 80706a1af2fSStefano Zampini bA->nnzstate[i*bA->nc+i] = subnnzstate; 80813135bc6SAlex Fikl } 80906a1af2fSStefano Zampini if (nnzstate) A->nonzerostate++; 81013135bc6SAlex Fikl PetscFunctionReturn(0); 81113135bc6SAlex Fikl } 81213135bc6SAlex Fikl 813f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx) 814f8170845SAlex Fikl { 815f8170845SAlex Fikl Mat_Nest *bA = (Mat_Nest*)A->data; 816f8170845SAlex Fikl PetscInt i,j; 817f8170845SAlex Fikl PetscErrorCode ierr; 818f8170845SAlex Fikl 819f8170845SAlex Fikl PetscFunctionBegin; 820f8170845SAlex Fikl for (i=0; i<bA->nr; i++) { 821f8170845SAlex Fikl for (j=0; j<bA->nc; j++) { 822f8170845SAlex Fikl if (bA->m[i][j]) { 823f8170845SAlex Fikl ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr); 824f8170845SAlex Fikl } 825f8170845SAlex Fikl } 826f8170845SAlex Fikl } 827f8170845SAlex Fikl PetscFunctionReturn(0); 828f8170845SAlex Fikl } 829f8170845SAlex Fikl 8302a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left) 831d8588912SDave May { 832d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 833d8588912SDave May Vec *L,*R; 834d8588912SDave May MPI_Comm comm; 835d8588912SDave May PetscInt i,j; 836d8588912SDave May PetscErrorCode ierr; 837d8588912SDave May 838d8588912SDave May PetscFunctionBegin; 839ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 840d8588912SDave May if (right) { 841d8588912SDave May /* allocate R */ 842854ce69bSBarry Smith ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr); 843d8588912SDave May /* Create the right vectors */ 844d8588912SDave May for (j=0; j<bA->nc; j++) { 845d8588912SDave May for (i=0; i<bA->nr; i++) { 846d8588912SDave May if (bA->m[i][j]) { 8472a7a6963SBarry Smith ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr); 848d8588912SDave May break; 849d8588912SDave May } 850d8588912SDave May } 8516c4ed002SBarry Smith if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column."); 852d8588912SDave May } 853f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr); 854d8588912SDave May /* hand back control to the nest vector */ 855d8588912SDave May for (j=0; j<bA->nc; j++) { 8566bf464f9SBarry Smith ierr = VecDestroy(&R[j]);CHKERRQ(ierr); 857d8588912SDave May } 858d8588912SDave May ierr = PetscFree(R);CHKERRQ(ierr); 859d8588912SDave May } 860d8588912SDave May 861d8588912SDave May if (left) { 862d8588912SDave May /* allocate L */ 863854ce69bSBarry Smith ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr); 864d8588912SDave May /* Create the left vectors */ 865d8588912SDave May for (i=0; i<bA->nr; i++) { 866d8588912SDave May for (j=0; j<bA->nc; j++) { 867d8588912SDave May if (bA->m[i][j]) { 8682a7a6963SBarry Smith ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr); 869d8588912SDave May break; 870d8588912SDave May } 871d8588912SDave May } 8726c4ed002SBarry Smith if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row."); 873d8588912SDave May } 874d8588912SDave May 875f349c1fdSJed Brown ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr); 876d8588912SDave May for (i=0; i<bA->nr; i++) { 8776bf464f9SBarry Smith ierr = VecDestroy(&L[i]);CHKERRQ(ierr); 878d8588912SDave May } 879d8588912SDave May 880d8588912SDave May ierr = PetscFree(L);CHKERRQ(ierr); 881d8588912SDave May } 882d8588912SDave May PetscFunctionReturn(0); 883d8588912SDave May } 884d8588912SDave May 885207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer) 886d8588912SDave May { 887d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 88829e60adbSStefano Zampini PetscBool isascii,viewSub = PETSC_FALSE; 889d8588912SDave May PetscInt i,j; 890d8588912SDave May PetscErrorCode ierr; 891d8588912SDave May 892d8588912SDave May PetscFunctionBegin; 893251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 894d8588912SDave May if (isascii) { 895d8588912SDave May 89629e60adbSStefano Zampini ierr = PetscOptionsGetBool(((PetscObject)A)->options,((PetscObject)A)->prefix,"-mat_view_nest_sub",&viewSub,NULL);CHKERRQ(ierr); 897d86155a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr); 898d86155a6SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 899d86155a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr); 900d8588912SDave May 901d86155a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr); 902d8588912SDave May for (i=0; i<bA->nr; i++) { 903d8588912SDave May for (j=0; j<bA->nc; j++) { 90419fd82e9SBarry Smith MatType type; 905270f95d7SJed Brown char name[256] = "",prefix[256] = ""; 906d8588912SDave May PetscInt NR,NC; 907d8588912SDave May PetscBool isNest = PETSC_FALSE; 908d8588912SDave May 909d8588912SDave May if (!bA->m[i][j]) { 910d86155a6SBarry Smith CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr); 911d8588912SDave May continue; 912d8588912SDave May } 913d8588912SDave May ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr); 914d8588912SDave May ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr); 9158caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);} 9168caf3d72SBarry Smith if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);} 917251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr); 918d8588912SDave May 919270f95d7SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr); 920d8588912SDave May 92129e60adbSStefano Zampini if (isNest || viewSub) { 922270f95d7SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* push1 */ 923d8588912SDave May ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr); 924270f95d7SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop1 */ 925d8588912SDave May } 926d8588912SDave May } 927d8588912SDave May } 928d86155a6SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); /* pop0 */ 929d8588912SDave May } 930d8588912SDave May PetscFunctionReturn(0); 931d8588912SDave May } 932d8588912SDave May 933207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A) 934d8588912SDave May { 935d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 936d8588912SDave May PetscInt i,j; 937d8588912SDave May PetscErrorCode ierr; 938d8588912SDave May 939d8588912SDave May PetscFunctionBegin; 940d8588912SDave May for (i=0; i<bA->nr; i++) { 941d8588912SDave May for (j=0; j<bA->nc; j++) { 942d8588912SDave May if (!bA->m[i][j]) continue; 943d8588912SDave May ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr); 944d8588912SDave May } 945d8588912SDave May } 946d8588912SDave May PetscFunctionReturn(0); 947d8588912SDave May } 948d8588912SDave May 949c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str) 950c222c20dSDavid Ham { 951c222c20dSDavid Ham Mat_Nest *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data; 952c222c20dSDavid Ham PetscInt i,j,nr = bA->nr,nc = bA->nc; 953c222c20dSDavid Ham PetscErrorCode ierr; 95406a1af2fSStefano Zampini PetscBool nnzstate = PETSC_FALSE; 955c222c20dSDavid Ham 956c222c20dSDavid Ham PetscFunctionBegin; 957c222c20dSDavid 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); 958c222c20dSDavid Ham for (i=0; i<nr; i++) { 959c222c20dSDavid Ham for (j=0; j<nc; j++) { 96006a1af2fSStefano Zampini PetscObjectState subnnzstate = 0; 96146a2b97cSJed Brown if (bA->m[i][j] && bB->m[i][j]) { 962c222c20dSDavid Ham ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr); 96346a2b97cSJed 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); 96406a1af2fSStefano Zampini ierr = MatGetNonzeroState(bB->m[i][j],&subnnzstate);CHKERRQ(ierr); 96506a1af2fSStefano Zampini nnzstate = (PetscBool)(nnzstate || bB->nnzstate[i*nc+j] != subnnzstate); 96606a1af2fSStefano Zampini bB->nnzstate[i*nc+j] = subnnzstate; 967c222c20dSDavid Ham } 968c222c20dSDavid Ham } 96906a1af2fSStefano Zampini if (nnzstate) B->nonzerostate++; 970c222c20dSDavid Ham PetscFunctionReturn(0); 971c222c20dSDavid Ham } 972c222c20dSDavid Ham 9736e76ffeaSPierre Jolivet static PetscErrorCode MatAXPY_Nest(Mat Y,PetscScalar a,Mat X,MatStructure str) 9746e76ffeaSPierre Jolivet { 9756e76ffeaSPierre Jolivet Mat_Nest *bY = (Mat_Nest*)Y->data,*bX = (Mat_Nest*)X->data; 9766e76ffeaSPierre Jolivet PetscInt i,j,nr = bY->nr,nc = bY->nc; 9776e76ffeaSPierre Jolivet PetscErrorCode ierr; 97806a1af2fSStefano Zampini PetscBool nnzstate = PETSC_FALSE; 9796e76ffeaSPierre Jolivet 9806e76ffeaSPierre Jolivet PetscFunctionBegin; 9816e76ffeaSPierre Jolivet if (nr != bX->nr || nc != bX->nc) SETERRQ4(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Cannot AXPY a MatNest of block size (%D,%D) with a MatNest of block size (%D,%D)",bX->nr,bX->nc,nr,nc); 9826e76ffeaSPierre Jolivet for (i=0; i<nr; i++) { 9836e76ffeaSPierre Jolivet for (j=0; j<nc; j++) { 98406a1af2fSStefano Zampini PetscObjectState subnnzstate = 0; 9856e76ffeaSPierre Jolivet if (bY->m[i][j] && bX->m[i][j]) { 9866e76ffeaSPierre Jolivet ierr = MatAXPY(bY->m[i][j],a,bX->m[i][j],str);CHKERRQ(ierr); 987c066aebcSStefano Zampini } else if (bX->m[i][j]) { 988c066aebcSStefano Zampini Mat M; 989c066aebcSStefano Zampini 990060bfc19SStefano Zampini if (str != DIFFERENT_NONZERO_PATTERN) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D. Use DIFFERENT_NONZERO_PATTERN",i,j); 991c066aebcSStefano Zampini ierr = MatDuplicate(bX->m[i][j],MAT_COPY_VALUES,&M);CHKERRQ(ierr); 992c066aebcSStefano Zampini ierr = MatNestSetSubMat(Y,i,j,M);CHKERRQ(ierr); 993c066aebcSStefano Zampini ierr = MatDestroy(&M);CHKERRQ(ierr); 994c066aebcSStefano Zampini } 995060bfc19SStefano Zampini if (bY->m[i][j]) { ierr = MatGetNonzeroState(bY->m[i][j],&subnnzstate);CHKERRQ(ierr); } 99606a1af2fSStefano Zampini nnzstate = (PetscBool)(nnzstate || bY->nnzstate[i*nc+j] != subnnzstate); 99706a1af2fSStefano Zampini bY->nnzstate[i*nc+j] = subnnzstate; 9986e76ffeaSPierre Jolivet } 9996e76ffeaSPierre Jolivet } 100006a1af2fSStefano Zampini if (nnzstate) Y->nonzerostate++; 10016e76ffeaSPierre Jolivet PetscFunctionReturn(0); 10026e76ffeaSPierre Jolivet } 10036e76ffeaSPierre Jolivet 1004207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B) 1005d8588912SDave May { 1006d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 1007841e96a3SJed Brown Mat *b; 1008841e96a3SJed Brown PetscInt i,j,nr = bA->nr,nc = bA->nc; 1009d8588912SDave May PetscErrorCode ierr; 1010d8588912SDave May 1011d8588912SDave May PetscFunctionBegin; 1012785e854fSJed Brown ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr); 1013841e96a3SJed Brown for (i=0; i<nr; i++) { 1014841e96a3SJed Brown for (j=0; j<nc; j++) { 1015841e96a3SJed Brown if (bA->m[i][j]) { 1016841e96a3SJed Brown ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr); 1017841e96a3SJed Brown } else { 10180298fd71SBarry Smith b[i*nc+j] = NULL; 1019d8588912SDave May } 1020d8588912SDave May } 1021d8588912SDave May } 1022ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr); 1023841e96a3SJed Brown /* Give the new MatNest exclusive ownership */ 1024841e96a3SJed Brown for (i=0; i<nr*nc; i++) { 10256bf464f9SBarry Smith ierr = MatDestroy(&b[i]);CHKERRQ(ierr); 1026d8588912SDave May } 1027d8588912SDave May ierr = PetscFree(b);CHKERRQ(ierr); 1028d8588912SDave May 1029841e96a3SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1030841e96a3SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1031d8588912SDave May PetscFunctionReturn(0); 1032d8588912SDave May } 1033d8588912SDave May 1034d8588912SDave May /* nest api */ 1035d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat) 1036d8588912SDave May { 1037d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 10385fd66863SKarl Rupp 1039d8588912SDave May PetscFunctionBegin; 1040ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 1041ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 1042d8588912SDave May *mat = bA->m[idxm][jdxm]; 1043d8588912SDave May PetscFunctionReturn(0); 1044d8588912SDave May } 1045d8588912SDave May 10469ba0d327SJed Brown /*@ 1047d8588912SDave May MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix. 1048d8588912SDave May 1049d8588912SDave May Not collective 1050d8588912SDave May 1051d8588912SDave May Input Parameters: 1052629881c0SJed Brown + A - nest matrix 1053d8588912SDave May . idxm - index of the matrix within the nest matrix 1054629881c0SJed Brown - jdxm - index of the matrix within the nest matrix 1055d8588912SDave May 1056d8588912SDave May Output Parameter: 1057d8588912SDave May . sub - matrix at index idxm,jdxm within the nest matrix 1058d8588912SDave May 1059d8588912SDave May Level: developer 1060d8588912SDave May 1061bb97c47cSPierre Jolivet .seealso: MatNestGetSize(), MatNestGetSubMats(), MatCreateNest(), MATNEST, MatNestSetSubMat(), 106279798668SBarry Smith MatNestGetLocalISs(), MatNestGetISs() 1063d8588912SDave May @*/ 10647087cfbeSBarry Smith PetscErrorCode MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub) 1065d8588912SDave May { 1066699a902aSJed Brown PetscErrorCode ierr; 1067d8588912SDave May 1068d8588912SDave May PetscFunctionBegin; 1069699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr); 1070d8588912SDave May PetscFunctionReturn(0); 1071d8588912SDave May } 1072d8588912SDave May 10730782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat) 10740782ca92SJed Brown { 10750782ca92SJed Brown Mat_Nest *bA = (Mat_Nest*)A->data; 10760782ca92SJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 10770782ca92SJed Brown PetscErrorCode ierr; 10780782ca92SJed Brown 10790782ca92SJed Brown PetscFunctionBegin; 1080ce94432eSBarry Smith if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1); 1081ce94432eSBarry Smith if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1); 10820782ca92SJed Brown ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 10830782ca92SJed Brown ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 10840782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr); 10850782ca92SJed Brown ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr); 10860782ca92SJed Brown ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr); 10870782ca92SJed Brown ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr); 1088ce94432eSBarry 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); 1089ce94432eSBarry 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); 109026fbe8dcSKarl Rupp 109106a1af2fSStefano Zampini /* do not increase object state */ 109206a1af2fSStefano Zampini if (mat == bA->m[idxm][jdxm]) PetscFunctionReturn(0); 109306a1af2fSStefano Zampini 10940782ca92SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 10950782ca92SJed Brown ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr); 10960782ca92SJed Brown bA->m[idxm][jdxm] = mat; 109706a1af2fSStefano Zampini ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 109806a1af2fSStefano Zampini ierr = MatGetNonzeroState(mat,&bA->nnzstate[idxm*bA->nc+jdxm]);CHKERRQ(ierr); 109906a1af2fSStefano Zampini A->nonzerostate++; 11000782ca92SJed Brown PetscFunctionReturn(0); 11010782ca92SJed Brown } 11020782ca92SJed Brown 11039ba0d327SJed Brown /*@ 11040782ca92SJed Brown MatNestSetSubMat - Set a single submatrix in the nest matrix. 11050782ca92SJed Brown 11060782ca92SJed Brown Logically collective on the submatrix communicator 11070782ca92SJed Brown 11080782ca92SJed Brown Input Parameters: 11090782ca92SJed Brown + A - nest matrix 11100782ca92SJed Brown . idxm - index of the matrix within the nest matrix 11110782ca92SJed Brown . jdxm - index of the matrix within the nest matrix 11120782ca92SJed Brown - sub - matrix at index idxm,jdxm within the nest matrix 11130782ca92SJed Brown 11140782ca92SJed Brown Notes: 11150782ca92SJed Brown The new submatrix must have the same size and communicator as that block of the nest. 11160782ca92SJed Brown 11170782ca92SJed Brown This increments the reference count of the submatrix. 11180782ca92SJed Brown 11190782ca92SJed Brown Level: developer 11200782ca92SJed Brown 1121bb97c47cSPierre Jolivet .seealso: MatNestSetSubMats(), MatNestGetSubMats(), MatNestGetLocalISs(), MATNEST, MatCreateNest(), 112279798668SBarry Smith MatNestGetSubMat(), MatNestGetISs(), MatNestGetSize() 11230782ca92SJed Brown @*/ 11240782ca92SJed Brown PetscErrorCode MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub) 11250782ca92SJed Brown { 11260782ca92SJed Brown PetscErrorCode ierr; 11270782ca92SJed Brown 11280782ca92SJed Brown PetscFunctionBegin; 11290782ca92SJed Brown ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr); 11300782ca92SJed Brown PetscFunctionReturn(0); 11310782ca92SJed Brown } 11320782ca92SJed Brown 1133d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 1134d8588912SDave May { 1135d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 11365fd66863SKarl Rupp 1137d8588912SDave May PetscFunctionBegin; 113826fbe8dcSKarl Rupp if (M) *M = bA->nr; 113926fbe8dcSKarl Rupp if (N) *N = bA->nc; 114026fbe8dcSKarl Rupp if (mat) *mat = bA->m; 1141d8588912SDave May PetscFunctionReturn(0); 1142d8588912SDave May } 1143d8588912SDave May 1144d8588912SDave May /*@C 1145d8588912SDave May MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix. 1146d8588912SDave May 1147d8588912SDave May Not collective 1148d8588912SDave May 1149d8588912SDave May Input Parameters: 1150629881c0SJed Brown . A - nest matrix 1151d8588912SDave May 1152d8588912SDave May Output Parameter: 1153629881c0SJed Brown + M - number of rows in the nest matrix 1154d8588912SDave May . N - number of cols in the nest matrix 1155629881c0SJed Brown - mat - 2d array of matrices 1156d8588912SDave May 1157d8588912SDave May Notes: 1158d8588912SDave May 1159d8588912SDave May The user should not free the array mat. 1160d8588912SDave May 1161351962e3SVincent Le Chenadec In Fortran, this routine has a calling sequence 1162351962e3SVincent Le Chenadec $ call MatNestGetSubMats(A, M, N, mat, ierr) 1163351962e3SVincent Le Chenadec where the space allocated for the optional argument mat is assumed large enough (if provided). 1164351962e3SVincent Le Chenadec 1165d8588912SDave May Level: developer 1166d8588912SDave May 1167bb97c47cSPierre Jolivet .seealso: MatNestGetSize(), MatNestGetSubMat(), MatNestGetLocalISs(), MATNEST, MatCreateNest(), 116879798668SBarry Smith MatNestSetSubMats(), MatNestGetISs(), MatNestSetSubMat() 1169d8588912SDave May @*/ 11707087cfbeSBarry Smith PetscErrorCode MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat) 1171d8588912SDave May { 1172699a902aSJed Brown PetscErrorCode ierr; 1173d8588912SDave May 1174d8588912SDave May PetscFunctionBegin; 1175699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr); 1176d8588912SDave May PetscFunctionReturn(0); 1177d8588912SDave May } 1178d8588912SDave May 11797087cfbeSBarry Smith PetscErrorCode MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N) 1180d8588912SDave May { 1181d8588912SDave May Mat_Nest *bA = (Mat_Nest*)A->data; 1182d8588912SDave May 1183d8588912SDave May PetscFunctionBegin; 118426fbe8dcSKarl Rupp if (M) *M = bA->nr; 118526fbe8dcSKarl Rupp if (N) *N = bA->nc; 1186d8588912SDave May PetscFunctionReturn(0); 1187d8588912SDave May } 1188d8588912SDave May 11899ba0d327SJed Brown /*@ 1190d8588912SDave May MatNestGetSize - Returns the size of the nest matrix. 1191d8588912SDave May 1192d8588912SDave May Not collective 1193d8588912SDave May 1194d8588912SDave May Input Parameters: 1195d8588912SDave May . A - nest matrix 1196d8588912SDave May 1197d8588912SDave May Output Parameter: 1198629881c0SJed Brown + M - number of rows in the nested mat 1199629881c0SJed Brown - N - number of cols in the nested mat 1200d8588912SDave May 1201d8588912SDave May Notes: 1202d8588912SDave May 1203d8588912SDave May Level: developer 1204d8588912SDave May 1205bb97c47cSPierre Jolivet .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MATNEST, MatCreateNest(), MatNestGetLocalISs(), 120679798668SBarry Smith MatNestGetISs() 1207d8588912SDave May @*/ 12087087cfbeSBarry Smith PetscErrorCode MatNestGetSize(Mat A,PetscInt *M,PetscInt *N) 1209d8588912SDave May { 1210699a902aSJed Brown PetscErrorCode ierr; 1211d8588912SDave May 1212d8588912SDave May PetscFunctionBegin; 1213699a902aSJed Brown ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr); 1214d8588912SDave May PetscFunctionReturn(0); 1215d8588912SDave May } 1216d8588912SDave May 1217f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[]) 1218900e7ff2SJed Brown { 1219900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 1220900e7ff2SJed Brown PetscInt i; 1221900e7ff2SJed Brown 1222900e7ff2SJed Brown PetscFunctionBegin; 1223900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i]; 1224900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i]; 1225900e7ff2SJed Brown PetscFunctionReturn(0); 1226900e7ff2SJed Brown } 1227900e7ff2SJed Brown 12283a4d7b9aSSatish Balay /*@C 1229900e7ff2SJed Brown MatNestGetISs - Returns the index sets partitioning the row and column spaces 1230900e7ff2SJed Brown 1231900e7ff2SJed Brown Not collective 1232900e7ff2SJed Brown 1233900e7ff2SJed Brown Input Parameters: 1234900e7ff2SJed Brown . A - nest matrix 1235900e7ff2SJed Brown 1236900e7ff2SJed Brown Output Parameter: 1237900e7ff2SJed Brown + rows - array of row index sets 1238900e7ff2SJed Brown - cols - array of column index sets 1239900e7ff2SJed Brown 1240900e7ff2SJed Brown Level: advanced 1241900e7ff2SJed Brown 1242900e7ff2SJed Brown Notes: 1243900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 1244900e7ff2SJed Brown 124579798668SBarry Smith .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs(), MATNEST, 1246bb97c47cSPierre Jolivet MatCreateNest(), MatNestGetSubMats(), MatNestSetSubMats() 1247900e7ff2SJed Brown @*/ 1248900e7ff2SJed Brown PetscErrorCode MatNestGetISs(Mat A,IS rows[],IS cols[]) 1249900e7ff2SJed Brown { 1250900e7ff2SJed Brown PetscErrorCode ierr; 1251900e7ff2SJed Brown 1252900e7ff2SJed Brown PetscFunctionBegin; 1253900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1254900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 1255900e7ff2SJed Brown PetscFunctionReturn(0); 1256900e7ff2SJed Brown } 1257900e7ff2SJed Brown 1258f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[]) 1259900e7ff2SJed Brown { 1260900e7ff2SJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 1261900e7ff2SJed Brown PetscInt i; 1262900e7ff2SJed Brown 1263900e7ff2SJed Brown PetscFunctionBegin; 1264900e7ff2SJed Brown if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i]; 1265900e7ff2SJed Brown if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i]; 1266900e7ff2SJed Brown PetscFunctionReturn(0); 1267900e7ff2SJed Brown } 1268900e7ff2SJed Brown 1269900e7ff2SJed Brown /*@C 1270900e7ff2SJed Brown MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces 1271900e7ff2SJed Brown 1272900e7ff2SJed Brown Not collective 1273900e7ff2SJed Brown 1274900e7ff2SJed Brown Input Parameters: 1275900e7ff2SJed Brown . A - nest matrix 1276900e7ff2SJed Brown 1277900e7ff2SJed Brown Output Parameter: 12780298fd71SBarry Smith + rows - array of row index sets (or NULL to ignore) 12790298fd71SBarry Smith - cols - array of column index sets (or NULL to ignore) 1280900e7ff2SJed Brown 1281900e7ff2SJed Brown Level: advanced 1282900e7ff2SJed Brown 1283900e7ff2SJed Brown Notes: 1284900e7ff2SJed Brown The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs. 1285900e7ff2SJed Brown 1286bb97c47cSPierre Jolivet .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs(), MatCreateNest(), 128779798668SBarry Smith MATNEST, MatNestSetSubMats(), MatNestSetSubMat() 1288900e7ff2SJed Brown @*/ 1289900e7ff2SJed Brown PetscErrorCode MatNestGetLocalISs(Mat A,IS rows[],IS cols[]) 1290900e7ff2SJed Brown { 1291900e7ff2SJed Brown PetscErrorCode ierr; 1292900e7ff2SJed Brown 1293900e7ff2SJed Brown PetscFunctionBegin; 1294900e7ff2SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1295900e7ff2SJed Brown ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr); 1296900e7ff2SJed Brown PetscFunctionReturn(0); 1297900e7ff2SJed Brown } 1298900e7ff2SJed Brown 129919fd82e9SBarry Smith PetscErrorCode MatNestSetVecType_Nest(Mat A,VecType vtype) 1300207556f9SJed Brown { 1301207556f9SJed Brown PetscErrorCode ierr; 1302207556f9SJed Brown PetscBool flg; 1303207556f9SJed Brown 1304207556f9SJed Brown PetscFunctionBegin; 1305207556f9SJed Brown ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr); 1306207556f9SJed Brown /* In reality, this only distinguishes VECNEST and "other" */ 13072a7a6963SBarry Smith if (flg) A->ops->getvecs = MatCreateVecs_Nest; 130812b53f24SSatish Balay else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0; 1309207556f9SJed Brown PetscFunctionReturn(0); 1310207556f9SJed Brown } 1311207556f9SJed Brown 1312207556f9SJed Brown /*@C 13132a7a6963SBarry Smith MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs() 1314207556f9SJed Brown 1315207556f9SJed Brown Not collective 1316207556f9SJed Brown 1317207556f9SJed Brown Input Parameters: 1318207556f9SJed Brown + A - nest matrix 1319207556f9SJed Brown - vtype - type to use for creating vectors 1320207556f9SJed Brown 1321207556f9SJed Brown Notes: 1322207556f9SJed Brown 1323207556f9SJed Brown Level: developer 1324207556f9SJed Brown 1325bb97c47cSPierre Jolivet .seealso: MatCreateVecs(), MATNEST, MatCreateNest() 1326207556f9SJed Brown @*/ 132719fd82e9SBarry Smith PetscErrorCode MatNestSetVecType(Mat A,VecType vtype) 1328207556f9SJed Brown { 1329207556f9SJed Brown PetscErrorCode ierr; 1330207556f9SJed Brown 1331207556f9SJed Brown PetscFunctionBegin; 133219fd82e9SBarry Smith ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr); 1333207556f9SJed Brown PetscFunctionReturn(0); 1334207556f9SJed Brown } 1335207556f9SJed Brown 1336c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1337d8588912SDave May { 1338c8883902SJed Brown Mat_Nest *s = (Mat_Nest*)A->data; 1339c8883902SJed Brown PetscInt i,j,m,n,M,N; 1340d8588912SDave May PetscErrorCode ierr; 134106a1af2fSStefano Zampini PetscBool cong; 1342d8588912SDave May 1343d8588912SDave May PetscFunctionBegin; 134406a1af2fSStefano Zampini ierr = MatReset_Nest(A);CHKERRQ(ierr); 134506a1af2fSStefano Zampini 1346c8883902SJed Brown s->nr = nr; 1347c8883902SJed Brown s->nc = nc; 1348d8588912SDave May 1349c8883902SJed Brown /* Create space for submatrices */ 1350854ce69bSBarry Smith ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr); 1351c8883902SJed Brown for (i=0; i<nr; i++) { 1352854ce69bSBarry Smith ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr); 1353d8588912SDave May } 1354c8883902SJed Brown for (i=0; i<nr; i++) { 1355c8883902SJed Brown for (j=0; j<nc; j++) { 1356c8883902SJed Brown s->m[i][j] = a[i*nc+j]; 1357c8883902SJed Brown if (a[i*nc+j]) { 1358c8883902SJed Brown ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr); 1359d8588912SDave May } 1360d8588912SDave May } 1361d8588912SDave May } 1362d8588912SDave May 13638188e55aSJed Brown ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr); 1364d8588912SDave May 1365854ce69bSBarry Smith ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr); 1366854ce69bSBarry Smith ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr); 1367c8883902SJed Brown for (i=0; i<nr; i++) s->row_len[i]=-1; 1368c8883902SJed Brown for (j=0; j<nc; j++) s->col_len[j]=-1; 1369d8588912SDave May 137006a1af2fSStefano Zampini ierr = PetscCalloc1(nr*nc,&s->nnzstate);CHKERRQ(ierr); 137106a1af2fSStefano Zampini for (i=0; i<nr; i++) { 137206a1af2fSStefano Zampini for (j=0; j<nc; j++) { 137306a1af2fSStefano Zampini if (s->m[i][j]) { 137406a1af2fSStefano Zampini ierr = MatGetNonzeroState(s->m[i][j],&s->nnzstate[i*nc+j]);CHKERRQ(ierr); 137506a1af2fSStefano Zampini } 137606a1af2fSStefano Zampini } 137706a1af2fSStefano Zampini } 137806a1af2fSStefano Zampini 13798188e55aSJed Brown ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr); 1380d8588912SDave May 1381c8883902SJed Brown ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr); 1382c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr); 1383c8883902SJed Brown ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr); 1384c8883902SJed Brown ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr); 1385c8883902SJed Brown 1386c8883902SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1387c8883902SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1388c8883902SJed Brown 138906a1af2fSStefano Zampini /* disable operations that are not supported for non-square matrices, 139006a1af2fSStefano Zampini or matrices for which is_row != is_col */ 139106a1af2fSStefano Zampini ierr = MatHasCongruentLayouts(A,&cong);CHKERRQ(ierr); 139206a1af2fSStefano Zampini if (cong && nr != nc) cong = PETSC_FALSE; 139306a1af2fSStefano Zampini if (cong) { 139406a1af2fSStefano Zampini for (i = 0; cong && i < nr; i++) { 1395320466b0SStefano Zampini ierr = ISEqualUnsorted(s->isglobal.row[i],s->isglobal.col[i],&cong);CHKERRQ(ierr); 139606a1af2fSStefano Zampini } 139706a1af2fSStefano Zampini } 139806a1af2fSStefano Zampini if (!cong) { 1399381b8e50SStefano Zampini A->ops->missingdiagonal = NULL; 140006a1af2fSStefano Zampini A->ops->getdiagonal = NULL; 140106a1af2fSStefano Zampini A->ops->shift = NULL; 140206a1af2fSStefano Zampini A->ops->diagonalset = NULL; 140306a1af2fSStefano Zampini } 140406a1af2fSStefano Zampini 14051795a4d1SJed Brown ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr); 140606a1af2fSStefano Zampini ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 140706a1af2fSStefano Zampini A->nonzerostate++; 1408d8588912SDave May PetscFunctionReturn(0); 1409d8588912SDave May } 1410d8588912SDave May 1411c8883902SJed Brown /*@ 1412c8883902SJed Brown MatNestSetSubMats - Sets the nested submatrices 1413c8883902SJed Brown 1414c8883902SJed Brown Collective on Mat 1415c8883902SJed Brown 1416c8883902SJed Brown Input Parameter: 1417ffd6319bSRichard Tran Mills + A - nested matrix 1418c8883902SJed Brown . nr - number of nested row blocks 14190298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1420c8883902SJed Brown . nc - number of nested column blocks 14210298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 14220298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1423c8883902SJed Brown 142406a1af2fSStefano Zampini Notes: this always resets any submatrix information previously set 142506a1af2fSStefano Zampini 1426c8883902SJed Brown Level: advanced 1427c8883902SJed Brown 142879798668SBarry Smith .seealso: MatCreateNest(), MATNEST, MatNestSetSubMat(), MatNestGetSubMat(), MatNestGetSubMats() 1429c8883902SJed Brown @*/ 1430c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[]) 1431c8883902SJed Brown { 1432c8883902SJed Brown PetscErrorCode ierr; 143306a1af2fSStefano Zampini PetscInt i; 1434c8883902SJed Brown 1435c8883902SJed Brown PetscFunctionBegin; 1436c8883902SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1437ce94432eSBarry Smith if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative"); 1438c8883902SJed Brown if (nr && is_row) { 1439c8883902SJed Brown PetscValidPointer(is_row,3); 1440c8883902SJed Brown for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3); 1441c8883902SJed Brown } 1442ce94432eSBarry Smith if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative"); 14431664e352SJed Brown if (nc && is_col) { 1444c8883902SJed Brown PetscValidPointer(is_col,5); 14459b30a8f6SBarry Smith for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5); 1446c8883902SJed Brown } 144706a1af2fSStefano Zampini if (nr*nc > 0) PetscValidPointer(a,6); 1448c8883902SJed 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); 1449c8883902SJed Brown PetscFunctionReturn(0); 1450c8883902SJed Brown } 1451d8588912SDave May 145245b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog) 145377019fcaSJed Brown { 145477019fcaSJed Brown PetscErrorCode ierr; 145577019fcaSJed Brown PetscBool flg; 145677019fcaSJed Brown PetscInt i,j,m,mi,*ix; 145777019fcaSJed Brown 145877019fcaSJed Brown PetscFunctionBegin; 1459aea6d515SStefano Zampini *ltog = NULL; 146077019fcaSJed Brown for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) { 146177019fcaSJed Brown if (islocal[i]) { 1462aea6d515SStefano Zampini ierr = ISGetLocalSize(islocal[i],&mi);CHKERRQ(ierr); 146377019fcaSJed Brown flg = PETSC_TRUE; /* We found a non-trivial entry */ 146477019fcaSJed Brown } else { 1465aea6d515SStefano Zampini ierr = ISGetLocalSize(isglobal[i],&mi);CHKERRQ(ierr); 146677019fcaSJed Brown } 146777019fcaSJed Brown m += mi; 146877019fcaSJed Brown } 1469aea6d515SStefano Zampini if (!flg) PetscFunctionReturn(0); 1470aea6d515SStefano Zampini 1471785e854fSJed Brown ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr); 1472165cd838SBarry Smith for (i=0,m=0; i<n; i++) { 14730298fd71SBarry Smith ISLocalToGlobalMapping smap = NULL; 1474e108cb99SStefano Zampini Mat sub = NULL; 1475f6d38dbbSStefano Zampini PetscSF sf; 1476f6d38dbbSStefano Zampini PetscLayout map; 1477aea6d515SStefano Zampini const PetscInt *ix2; 147877019fcaSJed Brown 1479165cd838SBarry Smith if (!colflg) { 148077019fcaSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 148177019fcaSJed Brown } else { 148277019fcaSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 148377019fcaSJed Brown } 1484191fd14bSBarry Smith if (sub) { 1485191fd14bSBarry Smith if (!colflg) { 1486191fd14bSBarry Smith ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr); 1487191fd14bSBarry Smith } else { 1488191fd14bSBarry Smith ierr = MatGetLocalToGlobalMapping(sub,NULL,&smap);CHKERRQ(ierr); 1489191fd14bSBarry Smith } 1490191fd14bSBarry Smith } 149177019fcaSJed Brown /* 149277019fcaSJed Brown Now we need to extract the monolithic global indices that correspond to the given split global indices. 149377019fcaSJed 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. 149477019fcaSJed Brown */ 1495aea6d515SStefano Zampini ierr = ISGetIndices(isglobal[i],&ix2);CHKERRQ(ierr); 1496aea6d515SStefano Zampini if (islocal[i]) { 1497aea6d515SStefano Zampini PetscInt *ilocal,*iremote; 1498aea6d515SStefano Zampini PetscInt mil,nleaves; 1499aea6d515SStefano Zampini 1500aea6d515SStefano Zampini ierr = ISGetLocalSize(islocal[i],&mi);CHKERRQ(ierr); 1501aea6d515SStefano Zampini if (!smap) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing local to global map"); 1502aea6d515SStefano Zampini for (j=0; j<mi; j++) ix[m+j] = j; 1503aea6d515SStefano Zampini ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr); 1504aea6d515SStefano Zampini 1505aea6d515SStefano Zampini /* PetscSFSetGraphLayout does not like negative indices */ 1506aea6d515SStefano Zampini ierr = PetscMalloc2(mi,&ilocal,mi,&iremote);CHKERRQ(ierr); 1507aea6d515SStefano Zampini for (j=0, nleaves = 0; j<mi; j++) { 1508aea6d515SStefano Zampini if (ix[m+j] < 0) continue; 1509aea6d515SStefano Zampini ilocal[nleaves] = j; 1510aea6d515SStefano Zampini iremote[nleaves] = ix[m+j]; 1511aea6d515SStefano Zampini nleaves++; 1512aea6d515SStefano Zampini } 1513aea6d515SStefano Zampini ierr = ISGetLocalSize(isglobal[i],&mil);CHKERRQ(ierr); 1514aea6d515SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 1515aea6d515SStefano Zampini ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)A),&map);CHKERRQ(ierr); 1516aea6d515SStefano Zampini ierr = PetscLayoutSetLocalSize(map,mil);CHKERRQ(ierr); 1517f6d38dbbSStefano Zampini ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1518aea6d515SStefano Zampini ierr = PetscSFSetGraphLayout(sf,map,nleaves,ilocal,PETSC_USE_POINTER,iremote);CHKERRQ(ierr); 1519f6d38dbbSStefano Zampini ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1520f6d38dbbSStefano Zampini ierr = PetscSFBcastBegin(sf,MPIU_INT,ix2,ix + m);CHKERRQ(ierr); 1521f6d38dbbSStefano Zampini ierr = PetscSFBcastEnd(sf,MPIU_INT,ix2,ix + m);CHKERRQ(ierr); 1522f6d38dbbSStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1523aea6d515SStefano Zampini ierr = PetscFree2(ilocal,iremote);CHKERRQ(ierr); 1524aea6d515SStefano Zampini } else { 1525aea6d515SStefano Zampini ierr = ISGetLocalSize(isglobal[i],&mi);CHKERRQ(ierr); 1526aea6d515SStefano Zampini for (j=0; j<mi; j++) ix[m+j] = ix2[i]; 1527aea6d515SStefano Zampini } 1528aea6d515SStefano Zampini ierr = ISRestoreIndices(isglobal[i],&ix2);CHKERRQ(ierr); 152977019fcaSJed Brown m += mi; 153077019fcaSJed Brown } 1531f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 153277019fcaSJed Brown PetscFunctionReturn(0); 153377019fcaSJed Brown } 153477019fcaSJed Brown 153577019fcaSJed Brown 1536d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */ 1537d8588912SDave May /* 1538d8588912SDave May nprocessors = NP 1539d8588912SDave May Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1)) 1540d8588912SDave May proc 0: => (g_0,h_0,) 1541d8588912SDave May proc 1: => (g_1,h_1,) 1542d8588912SDave May ... 1543d8588912SDave May proc nprocs-1: => (g_NP-1,h_NP-1,) 1544d8588912SDave May 1545d8588912SDave May proc 0: proc 1: proc nprocs-1: 1546d8588912SDave May is[0] = (0,1,2,...,nlocal(g_0)-1) (0,1,...,nlocal(g_1)-1) (0,1,...,nlocal(g_NP-1)) 1547d8588912SDave May 1548d8588912SDave May proc 0: 1549d8588912SDave May is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1) 1550d8588912SDave May proc 1: 1551d8588912SDave May is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1) 1552d8588912SDave May 1553d8588912SDave May proc NP-1: 1554d8588912SDave May is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1) 1555d8588912SDave May */ 1556841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[]) 1557d8588912SDave May { 1558e2d7f03fSJed Brown Mat_Nest *vs = (Mat_Nest*)A->data; 15598188e55aSJed Brown PetscInt i,j,offset,n,nsum,bs; 1560d8588912SDave May PetscErrorCode ierr; 15610298fd71SBarry Smith Mat sub = NULL; 1562d8588912SDave May 1563d8588912SDave May PetscFunctionBegin; 1564854ce69bSBarry Smith ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr); 1565854ce69bSBarry Smith ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr); 1566d8588912SDave May if (is_row) { /* valid IS is passed in */ 1567d8588912SDave May /* refs on is[] are incremeneted */ 1568e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1569d8588912SDave May ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr); 157026fbe8dcSKarl Rupp 1571e2d7f03fSJed Brown vs->isglobal.row[i] = is_row[i]; 1572d8588912SDave May } 15732ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each row */ 15748188e55aSJed Brown nsum = 0; 15758188e55aSJed Brown for (i=0; i<vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */ 15768188e55aSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 1577ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i); 15780298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 1579ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 15808188e55aSJed Brown nsum += n; 15818188e55aSJed Brown } 1582ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 158330bc264bSJed Brown offset -= nsum; 1584e2d7f03fSJed Brown for (i=0; i<vs->nr; i++) { 1585f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 15860298fd71SBarry Smith ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr); 15872ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1588ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr); 1589e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr); 15902ae74bdbSJed Brown offset += n; 1591d8588912SDave May } 1592d8588912SDave May } 1593d8588912SDave May 1594d8588912SDave May if (is_col) { /* valid IS is passed in */ 1595d8588912SDave May /* refs on is[] are incremeneted */ 1596e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1597d8588912SDave May ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr); 159826fbe8dcSKarl Rupp 1599e2d7f03fSJed Brown vs->isglobal.col[j] = is_col[j]; 1600d8588912SDave May } 16012ae74bdbSJed Brown } else { /* Create the ISs by inspecting sizes of a submatrix in each column */ 16022ae74bdbSJed Brown offset = A->cmap->rstart; 16038188e55aSJed Brown nsum = 0; 16048188e55aSJed Brown for (j=0; j<vs->nc; j++) { 16058188e55aSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 1606ce94432eSBarry Smith if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i); 16070298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 1608ce94432eSBarry Smith if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix"); 16098188e55aSJed Brown nsum += n; 16108188e55aSJed Brown } 1611ce94432eSBarry Smith ierr = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 161230bc264bSJed Brown offset -= nsum; 1613e2d7f03fSJed Brown for (j=0; j<vs->nc; j++) { 1614f349c1fdSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr); 16150298fd71SBarry Smith ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr); 16162ae74bdbSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1617ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr); 1618e2d7f03fSJed Brown ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr); 16192ae74bdbSJed Brown offset += n; 1620d8588912SDave May } 1621d8588912SDave May } 1622e2d7f03fSJed Brown 1623e2d7f03fSJed Brown /* Set up the local ISs */ 1624785e854fSJed Brown ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr); 1625785e854fSJed Brown ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr); 1626e2d7f03fSJed Brown for (i=0,offset=0; i<vs->nr; i++) { 1627e2d7f03fSJed Brown IS isloc; 16280298fd71SBarry Smith ISLocalToGlobalMapping rmap = NULL; 1629e2d7f03fSJed Brown PetscInt nlocal,bs; 1630e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr); 16310298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);} 1632207556f9SJed Brown if (rmap) { 1633e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1634e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr); 1635e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1636e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1637207556f9SJed Brown } else { 1638207556f9SJed Brown nlocal = 0; 16390298fd71SBarry Smith isloc = NULL; 1640207556f9SJed Brown } 1641e2d7f03fSJed Brown vs->islocal.row[i] = isloc; 1642e2d7f03fSJed Brown offset += nlocal; 1643e2d7f03fSJed Brown } 16448188e55aSJed Brown for (i=0,offset=0; i<vs->nc; i++) { 1645e2d7f03fSJed Brown IS isloc; 16460298fd71SBarry Smith ISLocalToGlobalMapping cmap = NULL; 1647e2d7f03fSJed Brown PetscInt nlocal,bs; 1648e2d7f03fSJed Brown ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr); 16490298fd71SBarry Smith if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);} 1650207556f9SJed Brown if (cmap) { 1651e2d7f03fSJed Brown ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr); 1652e2d7f03fSJed Brown ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr); 1653e2d7f03fSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr); 1654e2d7f03fSJed Brown ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr); 1655207556f9SJed Brown } else { 1656207556f9SJed Brown nlocal = 0; 16570298fd71SBarry Smith isloc = NULL; 1658207556f9SJed Brown } 1659e2d7f03fSJed Brown vs->islocal.col[i] = isloc; 1660e2d7f03fSJed Brown offset += nlocal; 1661e2d7f03fSJed Brown } 16620189643fSJed Brown 166377019fcaSJed Brown /* Set up the aggregate ISLocalToGlobalMapping */ 166477019fcaSJed Brown { 166545b6f7e9SBarry Smith ISLocalToGlobalMapping rmap,cmap; 166645b6f7e9SBarry Smith ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr); 166745b6f7e9SBarry Smith ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr); 166877019fcaSJed Brown if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);} 166977019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 167077019fcaSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 167177019fcaSJed Brown } 167277019fcaSJed Brown 16730189643fSJed Brown #if defined(PETSC_USE_DEBUG) 16740189643fSJed Brown for (i=0; i<vs->nr; i++) { 16750189643fSJed Brown for (j=0; j<vs->nc; j++) { 16760189643fSJed Brown PetscInt m,n,M,N,mi,ni,Mi,Ni; 16770189643fSJed Brown Mat B = vs->m[i][j]; 16780189643fSJed Brown if (!B) continue; 16790189643fSJed Brown ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr); 16800189643fSJed Brown ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr); 16810189643fSJed Brown ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr); 16820189643fSJed Brown ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr); 16830189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr); 16840189643fSJed Brown ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr); 1685ce94432eSBarry 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); 1686ce94432eSBarry 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); 16870189643fSJed Brown } 16880189643fSJed Brown } 16890189643fSJed Brown #endif 1690a061e289SJed Brown 1691a061e289SJed Brown /* Set A->assembled if all non-null blocks are currently assembled */ 1692a061e289SJed Brown for (i=0; i<vs->nr; i++) { 1693a061e289SJed Brown for (j=0; j<vs->nc; j++) { 1694a061e289SJed Brown if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0); 1695a061e289SJed Brown } 1696a061e289SJed Brown } 1697a061e289SJed Brown A->assembled = PETSC_TRUE; 1698d8588912SDave May PetscFunctionReturn(0); 1699d8588912SDave May } 1700d8588912SDave May 170145c38901SJed Brown /*@C 1702659c6bb0SJed Brown MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately 1703659c6bb0SJed Brown 1704659c6bb0SJed Brown Collective on Mat 1705659c6bb0SJed Brown 1706659c6bb0SJed Brown Input Parameter: 1707659c6bb0SJed Brown + comm - Communicator for the new Mat 1708659c6bb0SJed Brown . nr - number of nested row blocks 17090298fd71SBarry Smith . is_row - index sets for each nested row block, or NULL to make contiguous 1710659c6bb0SJed Brown . nc - number of nested column blocks 17110298fd71SBarry Smith . is_col - index sets for each nested column block, or NULL to make contiguous 17120298fd71SBarry Smith - a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL 1713659c6bb0SJed Brown 1714659c6bb0SJed Brown Output Parameter: 1715659c6bb0SJed Brown . B - new matrix 1716659c6bb0SJed Brown 1717659c6bb0SJed Brown Level: advanced 1718659c6bb0SJed Brown 171979798668SBarry Smith .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST, MatNestSetSubMat(), 172079798668SBarry Smith MatNestGetSubMat(), MatNestGetLocalISs(), MatNestGetSize(), 172179798668SBarry Smith MatNestGetISs(), MatNestSetSubMats(), MatNestGetSubMats() 1722659c6bb0SJed Brown @*/ 17237087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B) 1724d8588912SDave May { 1725d8588912SDave May Mat A; 1726d8588912SDave May PetscErrorCode ierr; 1727d8588912SDave May 1728d8588912SDave May PetscFunctionBegin; 1729c8883902SJed Brown *B = 0; 1730d8588912SDave May ierr = MatCreate(comm,&A);CHKERRQ(ierr); 1731c8883902SJed Brown ierr = MatSetType(A,MATNEST);CHKERRQ(ierr); 173291a28eb3SBarry Smith A->preallocated = PETSC_TRUE; 1733c8883902SJed Brown ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr); 1734d8588912SDave May *B = A; 1735d8588912SDave May PetscFunctionReturn(0); 1736d8588912SDave May } 1737659c6bb0SJed Brown 1738b68353e5Sstefano_zampini static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1739b68353e5Sstefano_zampini { 1740b68353e5Sstefano_zampini Mat_Nest *nest = (Mat_Nest*)A->data; 174123875855Sstefano_zampini Mat *trans; 1742b68353e5Sstefano_zampini PetscScalar **avv; 1743b68353e5Sstefano_zampini PetscScalar *vv; 1744b68353e5Sstefano_zampini PetscInt **aii,**ajj; 1745b68353e5Sstefano_zampini PetscInt *ii,*jj,*ci; 1746b68353e5Sstefano_zampini PetscInt nr,nc,nnz,i,j; 1747b68353e5Sstefano_zampini PetscBool done; 1748b68353e5Sstefano_zampini PetscErrorCode ierr; 1749b68353e5Sstefano_zampini 1750b68353e5Sstefano_zampini PetscFunctionBegin; 1751b68353e5Sstefano_zampini ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr); 1752b68353e5Sstefano_zampini if (reuse == MAT_REUSE_MATRIX) { 1753b68353e5Sstefano_zampini PetscInt rnr; 1754b68353e5Sstefano_zampini 1755b68353e5Sstefano_zampini ierr = MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);CHKERRQ(ierr); 1756b68353e5Sstefano_zampini if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ"); 1757b68353e5Sstefano_zampini if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows"); 1758b68353e5Sstefano_zampini ierr = MatSeqAIJGetArray(*newmat,&vv);CHKERRQ(ierr); 1759b68353e5Sstefano_zampini } 1760b68353e5Sstefano_zampini /* extract CSR for nested SeqAIJ matrices */ 1761b68353e5Sstefano_zampini nnz = 0; 176223875855Sstefano_zampini ierr = PetscCalloc4(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv,nest->nr*nest->nc,&trans);CHKERRQ(ierr); 1763b68353e5Sstefano_zampini for (i=0; i<nest->nr; ++i) { 1764b68353e5Sstefano_zampini for (j=0; j<nest->nc; ++j) { 1765b68353e5Sstefano_zampini Mat B = nest->m[i][j]; 1766b68353e5Sstefano_zampini if (B) { 1767b68353e5Sstefano_zampini PetscScalar *naa; 1768b68353e5Sstefano_zampini PetscInt *nii,*njj,nnr; 176923875855Sstefano_zampini PetscBool istrans; 1770b68353e5Sstefano_zampini 177123875855Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr); 177223875855Sstefano_zampini if (istrans) { 177323875855Sstefano_zampini Mat Bt; 177423875855Sstefano_zampini 177523875855Sstefano_zampini ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr); 177623875855Sstefano_zampini ierr = MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);CHKERRQ(ierr); 177723875855Sstefano_zampini B = trans[i*nest->nc+j]; 177823875855Sstefano_zampini } 1779b68353e5Sstefano_zampini ierr = MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);CHKERRQ(ierr); 1780b68353e5Sstefano_zampini if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ"); 1781b68353e5Sstefano_zampini ierr = MatSeqAIJGetArray(B,&naa);CHKERRQ(ierr); 1782b68353e5Sstefano_zampini nnz += nii[nnr]; 1783b68353e5Sstefano_zampini 1784b68353e5Sstefano_zampini aii[i*nest->nc+j] = nii; 1785b68353e5Sstefano_zampini ajj[i*nest->nc+j] = njj; 1786b68353e5Sstefano_zampini avv[i*nest->nc+j] = naa; 1787b68353e5Sstefano_zampini } 1788b68353e5Sstefano_zampini } 1789b68353e5Sstefano_zampini } 1790b68353e5Sstefano_zampini if (reuse != MAT_REUSE_MATRIX) { 1791b68353e5Sstefano_zampini ierr = PetscMalloc1(nr+1,&ii);CHKERRQ(ierr); 1792b68353e5Sstefano_zampini ierr = PetscMalloc1(nnz,&jj);CHKERRQ(ierr); 1793b68353e5Sstefano_zampini ierr = PetscMalloc1(nnz,&vv);CHKERRQ(ierr); 1794b68353e5Sstefano_zampini } else { 1795b68353e5Sstefano_zampini if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros"); 1796b68353e5Sstefano_zampini } 1797b68353e5Sstefano_zampini 1798b68353e5Sstefano_zampini /* new row pointer */ 1799580bdb30SBarry Smith ierr = PetscArrayzero(ii,nr+1);CHKERRQ(ierr); 1800b68353e5Sstefano_zampini for (i=0; i<nest->nr; ++i) { 1801b68353e5Sstefano_zampini PetscInt ncr,rst; 1802b68353e5Sstefano_zampini 1803b68353e5Sstefano_zampini ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr); 1804b68353e5Sstefano_zampini ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr); 1805b68353e5Sstefano_zampini for (j=0; j<nest->nc; ++j) { 1806b68353e5Sstefano_zampini if (aii[i*nest->nc+j]) { 1807b68353e5Sstefano_zampini PetscInt *nii = aii[i*nest->nc+j]; 1808b68353e5Sstefano_zampini PetscInt ir; 1809b68353e5Sstefano_zampini 1810b68353e5Sstefano_zampini for (ir=rst; ir<ncr+rst; ++ir) { 1811b68353e5Sstefano_zampini ii[ir+1] += nii[1]-nii[0]; 1812b68353e5Sstefano_zampini nii++; 1813b68353e5Sstefano_zampini } 1814b68353e5Sstefano_zampini } 1815b68353e5Sstefano_zampini } 1816b68353e5Sstefano_zampini } 1817b68353e5Sstefano_zampini for (i=0; i<nr; i++) ii[i+1] += ii[i]; 1818b68353e5Sstefano_zampini 1819b68353e5Sstefano_zampini /* construct CSR for the new matrix */ 1820b68353e5Sstefano_zampini ierr = PetscCalloc1(nr,&ci);CHKERRQ(ierr); 1821b68353e5Sstefano_zampini for (i=0; i<nest->nr; ++i) { 1822b68353e5Sstefano_zampini PetscInt ncr,rst; 1823b68353e5Sstefano_zampini 1824b68353e5Sstefano_zampini ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr); 1825b68353e5Sstefano_zampini ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr); 1826b68353e5Sstefano_zampini for (j=0; j<nest->nc; ++j) { 1827b68353e5Sstefano_zampini if (aii[i*nest->nc+j]) { 1828b68353e5Sstefano_zampini PetscScalar *nvv = avv[i*nest->nc+j]; 1829b68353e5Sstefano_zampini PetscInt *nii = aii[i*nest->nc+j]; 1830b68353e5Sstefano_zampini PetscInt *njj = ajj[i*nest->nc+j]; 1831b68353e5Sstefano_zampini PetscInt ir,cst; 1832b68353e5Sstefano_zampini 1833b68353e5Sstefano_zampini ierr = ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);CHKERRQ(ierr); 1834b68353e5Sstefano_zampini for (ir=rst; ir<ncr+rst; ++ir) { 1835b68353e5Sstefano_zampini PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir]; 1836b68353e5Sstefano_zampini 1837b68353e5Sstefano_zampini for (ij=0;ij<rsize;ij++) { 1838b68353e5Sstefano_zampini jj[ist+ij] = *njj+cst; 1839b68353e5Sstefano_zampini vv[ist+ij] = *nvv; 1840b68353e5Sstefano_zampini njj++; 1841b68353e5Sstefano_zampini nvv++; 1842b68353e5Sstefano_zampini } 1843b68353e5Sstefano_zampini ci[ir] += rsize; 1844b68353e5Sstefano_zampini nii++; 1845b68353e5Sstefano_zampini } 1846b68353e5Sstefano_zampini } 1847b68353e5Sstefano_zampini } 1848b68353e5Sstefano_zampini } 1849b68353e5Sstefano_zampini ierr = PetscFree(ci);CHKERRQ(ierr); 1850b68353e5Sstefano_zampini 1851b68353e5Sstefano_zampini /* restore info */ 1852b68353e5Sstefano_zampini for (i=0; i<nest->nr; ++i) { 1853b68353e5Sstefano_zampini for (j=0; j<nest->nc; ++j) { 1854b68353e5Sstefano_zampini Mat B = nest->m[i][j]; 1855b68353e5Sstefano_zampini if (B) { 1856b68353e5Sstefano_zampini PetscInt nnr = 0, k = i*nest->nc+j; 185723875855Sstefano_zampini 185823875855Sstefano_zampini B = (trans[k] ? trans[k] : B); 1859b68353e5Sstefano_zampini ierr = MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);CHKERRQ(ierr); 1860b68353e5Sstefano_zampini if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ"); 1861b68353e5Sstefano_zampini ierr = MatSeqAIJRestoreArray(B,&avv[k]);CHKERRQ(ierr); 186223875855Sstefano_zampini ierr = MatDestroy(&trans[k]);CHKERRQ(ierr); 1863b68353e5Sstefano_zampini } 1864b68353e5Sstefano_zampini } 1865b68353e5Sstefano_zampini } 186623875855Sstefano_zampini ierr = PetscFree4(aii,ajj,avv,trans);CHKERRQ(ierr); 1867b68353e5Sstefano_zampini 1868b68353e5Sstefano_zampini /* finalize newmat */ 1869b68353e5Sstefano_zampini if (reuse == MAT_INITIAL_MATRIX) { 1870b68353e5Sstefano_zampini ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);CHKERRQ(ierr); 1871b68353e5Sstefano_zampini } else if (reuse == MAT_INPLACE_MATRIX) { 1872b68353e5Sstefano_zampini Mat B; 1873b68353e5Sstefano_zampini 1874b68353e5Sstefano_zampini ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);CHKERRQ(ierr); 1875b68353e5Sstefano_zampini ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 1876b68353e5Sstefano_zampini } 1877b68353e5Sstefano_zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1878b68353e5Sstefano_zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1879b68353e5Sstefano_zampini { 1880b68353e5Sstefano_zampini Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data); 1881b68353e5Sstefano_zampini a->free_a = PETSC_TRUE; 1882b68353e5Sstefano_zampini a->free_ij = PETSC_TRUE; 1883b68353e5Sstefano_zampini } 1884b68353e5Sstefano_zampini PetscFunctionReturn(0); 1885b68353e5Sstefano_zampini } 1886b68353e5Sstefano_zampini 1887cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1888629c3df2SDmitry Karpeev { 1889629c3df2SDmitry Karpeev PetscErrorCode ierr; 1890629c3df2SDmitry Karpeev Mat_Nest *nest = (Mat_Nest*)A->data; 189183b1a929SMark Adams PetscInt m,n,M,N,i,j,k,*dnnz,*onnz,rstart; 1892649b366bSFande Kong PetscInt cstart,cend; 1893b68353e5Sstefano_zampini PetscMPIInt size; 1894629c3df2SDmitry Karpeev Mat C; 1895629c3df2SDmitry Karpeev 1896629c3df2SDmitry Karpeev PetscFunctionBegin; 1897b68353e5Sstefano_zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 1898b68353e5Sstefano_zampini if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */ 1899b68353e5Sstefano_zampini PetscInt nf; 1900b68353e5Sstefano_zampini PetscBool fast; 1901b68353e5Sstefano_zampini 1902b68353e5Sstefano_zampini ierr = PetscStrcmp(newtype,MATAIJ,&fast);CHKERRQ(ierr); 1903b68353e5Sstefano_zampini if (!fast) { 1904b68353e5Sstefano_zampini ierr = PetscStrcmp(newtype,MATSEQAIJ,&fast);CHKERRQ(ierr); 1905b68353e5Sstefano_zampini } 1906b68353e5Sstefano_zampini for (i=0; i<nest->nr && fast; ++i) { 1907b68353e5Sstefano_zampini for (j=0; j<nest->nc && fast; ++j) { 1908b68353e5Sstefano_zampini Mat B = nest->m[i][j]; 1909b68353e5Sstefano_zampini if (B) { 1910b68353e5Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);CHKERRQ(ierr); 191123875855Sstefano_zampini if (!fast) { 191223875855Sstefano_zampini PetscBool istrans; 191323875855Sstefano_zampini 191423875855Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr); 191523875855Sstefano_zampini if (istrans) { 191623875855Sstefano_zampini Mat Bt; 191723875855Sstefano_zampini 191823875855Sstefano_zampini ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr); 191923875855Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);CHKERRQ(ierr); 192023875855Sstefano_zampini } 1921b68353e5Sstefano_zampini } 1922b68353e5Sstefano_zampini } 1923b68353e5Sstefano_zampini } 1924b68353e5Sstefano_zampini } 1925b68353e5Sstefano_zampini for (i=0, nf=0; i<nest->nr && fast; ++i) { 1926b68353e5Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);CHKERRQ(ierr); 1927b68353e5Sstefano_zampini if (fast) { 1928b68353e5Sstefano_zampini PetscInt f,s; 1929b68353e5Sstefano_zampini 1930b68353e5Sstefano_zampini ierr = ISStrideGetInfo(nest->isglobal.row[i],&f,&s);CHKERRQ(ierr); 1931b68353e5Sstefano_zampini if (f != nf || s != 1) { fast = PETSC_FALSE; } 1932b68353e5Sstefano_zampini else { 1933b68353e5Sstefano_zampini ierr = ISGetSize(nest->isglobal.row[i],&f);CHKERRQ(ierr); 1934b68353e5Sstefano_zampini nf += f; 1935b68353e5Sstefano_zampini } 1936b68353e5Sstefano_zampini } 1937b68353e5Sstefano_zampini } 1938b68353e5Sstefano_zampini for (i=0, nf=0; i<nest->nc && fast; ++i) { 1939b68353e5Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);CHKERRQ(ierr); 1940b68353e5Sstefano_zampini if (fast) { 1941b68353e5Sstefano_zampini PetscInt f,s; 1942b68353e5Sstefano_zampini 1943b68353e5Sstefano_zampini ierr = ISStrideGetInfo(nest->isglobal.col[i],&f,&s);CHKERRQ(ierr); 1944b68353e5Sstefano_zampini if (f != nf || s != 1) { fast = PETSC_FALSE; } 1945b68353e5Sstefano_zampini else { 1946b68353e5Sstefano_zampini ierr = ISGetSize(nest->isglobal.col[i],&f);CHKERRQ(ierr); 1947b68353e5Sstefano_zampini nf += f; 1948b68353e5Sstefano_zampini } 1949b68353e5Sstefano_zampini } 1950b68353e5Sstefano_zampini } 1951b68353e5Sstefano_zampini if (fast) { 1952b68353e5Sstefano_zampini ierr = MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);CHKERRQ(ierr); 1953b68353e5Sstefano_zampini PetscFunctionReturn(0); 1954b68353e5Sstefano_zampini } 1955b68353e5Sstefano_zampini } 1956629c3df2SDmitry Karpeev ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1957629c3df2SDmitry Karpeev ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1958649b366bSFande Kong ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); 1959629c3df2SDmitry Karpeev switch (reuse) { 1960629c3df2SDmitry Karpeev case MAT_INITIAL_MATRIX: 1961ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 1962629c3df2SDmitry Karpeev ierr = MatSetType(C,newtype);CHKERRQ(ierr); 1963629c3df2SDmitry Karpeev ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 1964629c3df2SDmitry Karpeev *newmat = C; 1965629c3df2SDmitry Karpeev break; 1966629c3df2SDmitry Karpeev case MAT_REUSE_MATRIX: 1967629c3df2SDmitry Karpeev C = *newmat; 1968629c3df2SDmitry Karpeev break; 1969ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse"); 1970629c3df2SDmitry Karpeev } 1971785e854fSJed Brown ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr); 1972629c3df2SDmitry Karpeev onnz = dnnz + m; 1973629c3df2SDmitry Karpeev for (k=0; k<m; k++) { 1974629c3df2SDmitry Karpeev dnnz[k] = 0; 1975629c3df2SDmitry Karpeev onnz[k] = 0; 1976629c3df2SDmitry Karpeev } 1977629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 1978629c3df2SDmitry Karpeev IS bNis; 1979629c3df2SDmitry Karpeev PetscInt bN; 1980629c3df2SDmitry Karpeev const PetscInt *bNindices; 1981629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 1982629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 1983629c3df2SDmitry Karpeev ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr); 1984629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 1985629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 1986629c3df2SDmitry Karpeev PetscSF bmsf; 1987649b366bSFande Kong PetscSFNode *iremote; 1988629c3df2SDmitry Karpeev Mat B; 1989649b366bSFande Kong PetscInt bm, *sub_dnnz,*sub_onnz, br; 1990629c3df2SDmitry Karpeev const PetscInt *bmindices; 1991629c3df2SDmitry Karpeev B = nest->m[i][j]; 1992629c3df2SDmitry Karpeev if (!B) continue; 1993629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 1994629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 1995ce94432eSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr); 1996649b366bSFande Kong ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr); 1997649b366bSFande Kong ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr); 1998649b366bSFande Kong ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr); 1999649b366bSFande Kong for (k = 0; k < bm; ++k){ 2000649b366bSFande Kong sub_dnnz[k] = 0; 2001649b366bSFande Kong sub_onnz[k] = 0; 2002649b366bSFande Kong } 2003629c3df2SDmitry Karpeev /* 2004629c3df2SDmitry Karpeev Locate the owners for all of the locally-owned global row indices for this row block. 2005629c3df2SDmitry Karpeev These determine the roots of PetscSF used to communicate preallocation data to row owners. 2006629c3df2SDmitry Karpeev The roots correspond to the dnnz and onnz entries; thus, there are two roots per row. 2007629c3df2SDmitry Karpeev */ 200883b1a929SMark Adams ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 2009629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 2010131c27b5Sprj- PetscInt row = bmindices[br], brncols, col; 2011629c3df2SDmitry Karpeev const PetscInt *brcols; 2012a4b3d3acSMatthew G Knepley PetscInt rowrel = 0; /* row's relative index on its owner rank */ 2013131c27b5Sprj- PetscMPIInt rowowner = 0; 2014629c3df2SDmitry Karpeev ierr = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr); 2015649b366bSFande Kong /* how many roots */ 2016649b366bSFande Kong iremote[br].rank = rowowner; iremote[br].index = rowrel; /* edge from bmdnnz to dnnz */ 2017649b366bSFande Kong /* get nonzero pattern */ 201883b1a929SMark Adams ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr); 2019629c3df2SDmitry Karpeev for (k=0; k<brncols; k++) { 2020629c3df2SDmitry Karpeev col = bNindices[brcols[k]]; 2021649b366bSFande Kong if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) { 2022649b366bSFande Kong sub_dnnz[br]++; 2023649b366bSFande Kong } else { 2024649b366bSFande Kong sub_onnz[br]++; 2025649b366bSFande Kong } 2026629c3df2SDmitry Karpeev } 202783b1a929SMark Adams ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr); 2028629c3df2SDmitry Karpeev } 2029629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 2030629c3df2SDmitry Karpeev /* bsf will have to take care of disposing of bedges. */ 2031649b366bSFande Kong ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 2032649b366bSFande Kong ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr); 2033649b366bSFande Kong ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr); 2034649b366bSFande Kong ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr); 2035649b366bSFande Kong ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr); 2036649b366bSFande Kong ierr = PetscFree(sub_dnnz);CHKERRQ(ierr); 2037649b366bSFande Kong ierr = PetscFree(sub_onnz);CHKERRQ(ierr); 2038629c3df2SDmitry Karpeev ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr); 2039629c3df2SDmitry Karpeev } 204022d28d08SBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 2041629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 204265a4a0a3Sstefano_zampini } 204365a4a0a3Sstefano_zampini /* Resize preallocation if overestimated */ 204465a4a0a3Sstefano_zampini for (i=0;i<m;i++) { 204565a4a0a3Sstefano_zampini dnnz[i] = PetscMin(dnnz[i],A->cmap->n); 204665a4a0a3Sstefano_zampini onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n); 2047629c3df2SDmitry Karpeev } 2048629c3df2SDmitry Karpeev ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr); 2049629c3df2SDmitry Karpeev ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr); 2050629c3df2SDmitry Karpeev ierr = PetscFree(dnnz);CHKERRQ(ierr); 2051629c3df2SDmitry Karpeev 2052629c3df2SDmitry Karpeev /* Fill by row */ 2053629c3df2SDmitry Karpeev for (j=0; j<nest->nc; ++j) { 2054629c3df2SDmitry Karpeev /* Using global column indices and ISAllGather() is not scalable. */ 2055629c3df2SDmitry Karpeev IS bNis; 2056629c3df2SDmitry Karpeev PetscInt bN; 2057629c3df2SDmitry Karpeev const PetscInt *bNindices; 2058629c3df2SDmitry Karpeev ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr); 2059629c3df2SDmitry Karpeev ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr); 2060629c3df2SDmitry Karpeev ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr); 2061629c3df2SDmitry Karpeev for (i=0; i<nest->nr; ++i) { 2062629c3df2SDmitry Karpeev Mat B; 2063629c3df2SDmitry Karpeev PetscInt bm, br; 2064629c3df2SDmitry Karpeev const PetscInt *bmindices; 2065629c3df2SDmitry Karpeev B = nest->m[i][j]; 2066629c3df2SDmitry Karpeev if (!B) continue; 2067629c3df2SDmitry Karpeev ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr); 2068629c3df2SDmitry Karpeev ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 206983b1a929SMark Adams ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 2070629c3df2SDmitry Karpeev for (br = 0; br < bm; ++br) { 2071629c3df2SDmitry Karpeev PetscInt row = bmindices[br], brncols, *cols; 2072629c3df2SDmitry Karpeev const PetscInt *brcols; 2073629c3df2SDmitry Karpeev const PetscScalar *brcoldata; 207483b1a929SMark Adams ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 2075785e854fSJed Brown ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr); 207626fbe8dcSKarl Rupp for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]]; 2077629c3df2SDmitry Karpeev /* 2078629c3df2SDmitry Karpeev Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match. 2079629c3df2SDmitry Karpeev Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES. 2080629c3df2SDmitry Karpeev */ 2081a2ea699eSBarry Smith ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr); 208283b1a929SMark Adams ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr); 2083629c3df2SDmitry Karpeev ierr = PetscFree(cols);CHKERRQ(ierr); 2084629c3df2SDmitry Karpeev } 2085629c3df2SDmitry Karpeev ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr); 2086629c3df2SDmitry Karpeev } 2087a2ea699eSBarry Smith ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr); 2088629c3df2SDmitry Karpeev ierr = ISDestroy(&bNis);CHKERRQ(ierr); 2089629c3df2SDmitry Karpeev } 2090629c3df2SDmitry Karpeev ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2091629c3df2SDmitry Karpeev ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2092629c3df2SDmitry Karpeev PetscFunctionReturn(0); 2093629c3df2SDmitry Karpeev } 2094629c3df2SDmitry Karpeev 20958b7d3b4bSBarry Smith PetscErrorCode MatHasOperation_Nest(Mat mat,MatOperation op,PetscBool *has) 20968b7d3b4bSBarry Smith { 20978b7d3b4bSBarry Smith Mat_Nest *bA = (Mat_Nest*)mat->data; 20988b7d3b4bSBarry Smith PetscInt i,j,nr = bA->nr,nc = bA->nc; 20998b7d3b4bSBarry Smith PetscBool flg; 210052c5f739Sprj- PetscErrorCode ierr; 210152c5f739Sprj- PetscFunctionBegin; 21028b7d3b4bSBarry Smith 210352c5f739Sprj- *has = PETSC_FALSE; 210452c5f739Sprj- if (op == MATOP_MULT_TRANSPOSE || op == MATOP_MAT_MULT) { 21058b7d3b4bSBarry Smith for (j=0; j<nc; j++) { 21068b7d3b4bSBarry Smith for (i=0; i<nr; i++) { 21078b7d3b4bSBarry Smith if (!bA->m[i][j]) continue; 210852c5f739Sprj- ierr = MatHasOperation(bA->m[i][j],op,&flg);CHKERRQ(ierr); 21098b7d3b4bSBarry Smith if (!flg) PetscFunctionReturn(0); 21108b7d3b4bSBarry Smith } 21118b7d3b4bSBarry Smith } 21128b7d3b4bSBarry Smith } 211352c5f739Sprj- if (((void**)mat->ops)[op] || (op == MATOP_MAT_MULT && flg)) *has = PETSC_TRUE; 21148b7d3b4bSBarry Smith PetscFunctionReturn(0); 21158b7d3b4bSBarry Smith } 21168b7d3b4bSBarry Smith 2117659c6bb0SJed Brown /*MC 2118659c6bb0SJed Brown MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately. 2119659c6bb0SJed Brown 2120659c6bb0SJed Brown Level: intermediate 2121659c6bb0SJed Brown 2122659c6bb0SJed Brown Notes: 2123659c6bb0SJed Brown This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices. 2124659c6bb0SJed Brown It allows the use of symmetric and block formats for parts of multi-physics simulations. 2125950540a4SJed Brown It is usually used with DMComposite and DMCreateMatrix() 2126659c6bb0SJed Brown 21278b7d3b4bSBarry Smith Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero 21288b7d3b4bSBarry Smith rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes 21298b7d3b4bSBarry Smith than the nest matrix. 21308b7d3b4bSBarry Smith 213179798668SBarry Smith .seealso: MatCreate(), MatType, MatCreateNest(), MatNestSetSubMat(), MatNestGetSubMat(), 213279798668SBarry Smith VecCreateNest(), DMCreateMatrix(), DMCOMPOSITE, MatNestSetVecType(), MatNestGetLocalISs(), 213379798668SBarry Smith MatNestGetISs(), MatNestSetSubMats(), MatNestGetSubMats() 2134659c6bb0SJed Brown M*/ 21358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A) 2136c8883902SJed Brown { 2137c8883902SJed Brown Mat_Nest *s; 2138c8883902SJed Brown PetscErrorCode ierr; 2139c8883902SJed Brown 2140c8883902SJed Brown PetscFunctionBegin; 2141b00a9115SJed Brown ierr = PetscNewLog(A,&s);CHKERRQ(ierr); 2142c8883902SJed Brown A->data = (void*)s; 2143e7c19651SJed Brown 2144e7c19651SJed Brown s->nr = -1; 2145e7c19651SJed Brown s->nc = -1; 21460298fd71SBarry Smith s->m = NULL; 2147e7c19651SJed Brown s->splitassembly = PETSC_FALSE; 2148c8883902SJed Brown 2149c8883902SJed Brown ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr); 215026fbe8dcSKarl Rupp 2151c8883902SJed Brown A->ops->mult = MatMult_Nest; 21529194d70fSJed Brown A->ops->multadd = MatMultAdd_Nest; 2153c8883902SJed Brown A->ops->multtranspose = MatMultTranspose_Nest; 21549194d70fSJed Brown A->ops->multtransposeadd = MatMultTransposeAdd_Nest; 2155f8170845SAlex Fikl A->ops->transpose = MatTranspose_Nest; 2156c8883902SJed Brown A->ops->assemblybegin = MatAssemblyBegin_Nest; 2157c8883902SJed Brown A->ops->assemblyend = MatAssemblyEnd_Nest; 2158c8883902SJed Brown A->ops->zeroentries = MatZeroEntries_Nest; 2159c222c20dSDavid Ham A->ops->copy = MatCopy_Nest; 21606e76ffeaSPierre Jolivet A->ops->axpy = MatAXPY_Nest; 2161c8883902SJed Brown A->ops->duplicate = MatDuplicate_Nest; 21627dae84e0SHong Zhang A->ops->createsubmatrix = MatCreateSubMatrix_Nest; 2163c8883902SJed Brown A->ops->destroy = MatDestroy_Nest; 2164c8883902SJed Brown A->ops->view = MatView_Nest; 2165c8883902SJed Brown A->ops->getvecs = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */ 2166c8883902SJed Brown A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_Nest; 2167c8883902SJed Brown A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest; 2168429bac76SJed Brown A->ops->getdiagonal = MatGetDiagonal_Nest; 2169429bac76SJed Brown A->ops->diagonalscale = MatDiagonalScale_Nest; 2170a061e289SJed Brown A->ops->scale = MatScale_Nest; 2171a061e289SJed Brown A->ops->shift = MatShift_Nest; 217213135bc6SAlex Fikl A->ops->diagonalset = MatDiagonalSet_Nest; 2173f8170845SAlex Fikl A->ops->setrandom = MatSetRandom_Nest; 21748b7d3b4bSBarry Smith A->ops->hasoperation = MatHasOperation_Nest; 2175381b8e50SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_Nest; 2176c8883902SJed Brown 2177c8883902SJed Brown A->spptr = 0; 2178c8883902SJed Brown A->assembled = PETSC_FALSE; 2179c8883902SJed Brown 2180c8883902SJed Brown /* expose Nest api's */ 2181bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C", MatNestGetSubMat_Nest);CHKERRQ(ierr); 2182bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C", MatNestSetSubMat_Nest);CHKERRQ(ierr); 2183bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C", MatNestGetSubMats_Nest);CHKERRQ(ierr); 2184bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C", MatNestGetSize_Nest);CHKERRQ(ierr); 2185bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C", MatNestGetISs_Nest);CHKERRQ(ierr); 2186bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr); 2187bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C", MatNestSetVecType_Nest);CHKERRQ(ierr); 2188bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C", MatNestSetSubMats_Nest);CHKERRQ(ierr); 21890899c546SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C", MatConvert_Nest_AIJ);CHKERRQ(ierr); 21900899c546SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C", MatConvert_Nest_AIJ);CHKERRQ(ierr); 219183b1a929SMark Adams ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C", MatConvert_Nest_AIJ);CHKERRQ(ierr); 21925e3038f0Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C", MatConvert_Nest_IS);CHKERRQ(ierr); 21934222ddf1SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_seqdense_C",MatProductSetFromOptions_Nest_Dense);CHKERRQ(ierr); 21944222ddf1SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_mpidense_C",MatProductSetFromOptions_Nest_Dense);CHKERRQ(ierr); 21954222ddf1SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_dense_C",MatProductSetFromOptions_Nest_Dense);CHKERRQ(ierr); 2196c8883902SJed Brown 2197c8883902SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr); 2198c8883902SJed Brown PetscFunctionReturn(0); 2199c8883902SJed Brown } 2200