1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith We provide: 9b4319ba4SBarry Smith MatMult() 10b4319ba4SBarry Smith 11b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 12b4319ba4SBarry Smith 13b4319ba4SBarry Smith */ 14b4319ba4SBarry Smith 15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 16b4319ba4SBarry Smith 17b4319ba4SBarry Smith #undef __FUNCT__ 1869796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 1969796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 2069796d55SStefano Zampini { 2169796d55SStefano Zampini PetscErrorCode ierr; 2269796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 2369796d55SStefano Zampini PetscBool local_sym; 2469796d55SStefano Zampini 2569796d55SStefano Zampini PetscFunctionBegin; 2669796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 2769796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 2869796d55SStefano Zampini PetscFunctionReturn(0); 2969796d55SStefano Zampini } 3069796d55SStefano Zampini 3169796d55SStefano Zampini #undef __FUNCT__ 3269796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 3369796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 3469796d55SStefano Zampini { 3569796d55SStefano Zampini PetscErrorCode ierr; 3669796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 3769796d55SStefano Zampini PetscBool local_sym; 3869796d55SStefano Zampini 3969796d55SStefano Zampini PetscFunctionBegin; 4069796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 4169796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 4269796d55SStefano Zampini PetscFunctionReturn(0); 4369796d55SStefano Zampini } 4469796d55SStefano Zampini 4569796d55SStefano Zampini #undef __FUNCT__ 46b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 47dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 48b4319ba4SBarry Smith { 49dfbe8321SBarry Smith PetscErrorCode ierr; 50b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 51b4319ba4SBarry Smith 52b4319ba4SBarry Smith PetscFunctionBegin; 536bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 546bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 556bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 566bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 576bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 58bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 59dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 60bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 61b4319ba4SBarry Smith PetscFunctionReturn(0); 62b4319ba4SBarry Smith } 63b4319ba4SBarry Smith 64b4319ba4SBarry Smith #undef __FUNCT__ 65b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 66dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 67b4319ba4SBarry Smith { 68dfbe8321SBarry Smith PetscErrorCode ierr; 69b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 70b4319ba4SBarry Smith PetscScalar zero = 0.0; 71b4319ba4SBarry Smith 72b4319ba4SBarry Smith PetscFunctionBegin; 73b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 74ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 75ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 76b4319ba4SBarry Smith 77b4319ba4SBarry Smith /* multiply the local matrix */ 78b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 79b4319ba4SBarry Smith 80b4319ba4SBarry Smith /* scatter product back into global memory */ 812dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 82ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 83ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 84b4319ba4SBarry Smith PetscFunctionReturn(0); 85b4319ba4SBarry Smith } 86b4319ba4SBarry Smith 87b4319ba4SBarry Smith #undef __FUNCT__ 882e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 892e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 902e74eeadSLisandro Dalcin { 91650997f4SStefano Zampini Vec temp_vec; 922e74eeadSLisandro Dalcin PetscErrorCode ierr; 932e74eeadSLisandro Dalcin 942e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 95650997f4SStefano Zampini if (v3 != v2) { 96650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 97650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 98650997f4SStefano Zampini } else { 99650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 100650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 101650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 102650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 103650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 104650997f4SStefano Zampini } 1052e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1062e74eeadSLisandro Dalcin } 1072e74eeadSLisandro Dalcin 1082e74eeadSLisandro Dalcin #undef __FUNCT__ 1092e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 1102e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 1112e74eeadSLisandro Dalcin { 1122e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1132e74eeadSLisandro Dalcin PetscErrorCode ierr; 1142e74eeadSLisandro Dalcin 1152e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 1162e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 117ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 118ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1192e74eeadSLisandro Dalcin 1202e74eeadSLisandro Dalcin /* multiply the local matrix */ 1212e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 1222e74eeadSLisandro Dalcin 1232e74eeadSLisandro Dalcin /* scatter product back into global vector */ 1242e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 125ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 126ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1272e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1282e74eeadSLisandro Dalcin } 1292e74eeadSLisandro Dalcin 1302e74eeadSLisandro Dalcin #undef __FUNCT__ 1312e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1322e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 1332e74eeadSLisandro Dalcin { 134650997f4SStefano Zampini Vec temp_vec; 1352e74eeadSLisandro Dalcin PetscErrorCode ierr; 1362e74eeadSLisandro Dalcin 1372e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 138650997f4SStefano Zampini if (v3 != v2) { 139650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 140650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 141650997f4SStefano Zampini } else { 142650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 143650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 144650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 145650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 146650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 147650997f4SStefano Zampini } 1482e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1492e74eeadSLisandro Dalcin } 1502e74eeadSLisandro Dalcin 1512e74eeadSLisandro Dalcin #undef __FUNCT__ 152b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 153dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 154b4319ba4SBarry Smith { 155b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 156dfbe8321SBarry Smith PetscErrorCode ierr; 157b4319ba4SBarry Smith PetscViewer sviewer; 158b4319ba4SBarry Smith 159b4319ba4SBarry Smith PetscFunctionBegin; 160b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 161b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 162b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 163b4319ba4SBarry Smith PetscFunctionReturn(0); 164b4319ba4SBarry Smith } 165b4319ba4SBarry Smith 166b4319ba4SBarry Smith #undef __FUNCT__ 167b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 168784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 169b4319ba4SBarry Smith { 170dfbe8321SBarry Smith PetscErrorCode ierr; 1714e4c7dbeSStefano Zampini PetscInt n,bs; 172b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 173b4319ba4SBarry Smith IS from,to; 174b4319ba4SBarry Smith Vec global; 175b4319ba4SBarry Smith 176b4319ba4SBarry Smith PetscFunctionBegin; 177784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 178ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 179*70cf5478SStefano Zampini if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 180*70cf5478SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 181*70cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 182*70cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 183*70cf5478SStefano Zampini ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 184*70cf5478SStefano Zampini } 185784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 1866bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 187784ac674SJed Brown is->mapping = rmapping; 188fa7f1dd8SStefano Zampini /* 189fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 190fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 191fa7f1dd8SStefano Zampini */ 192b4319ba4SBarry Smith 193b4319ba4SBarry Smith /* Create the local matrix A */ 194784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 1954e4c7dbeSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 196f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 197f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 1984e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 199ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 200ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 201b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 202b4319ba4SBarry Smith 203b4319ba4SBarry Smith /* Create the local work vectors */ 2044e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 2054e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 2064e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 207ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 208ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 2094e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 210b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 211b4319ba4SBarry Smith 212b4319ba4SBarry Smith /* setup the global to local scatter */ 213b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 214784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 2150298fd71SBarry Smith ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr); 216b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 2176bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 2186bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 2196bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 220b4319ba4SBarry Smith PetscFunctionReturn(0); 221b4319ba4SBarry Smith } 222b4319ba4SBarry Smith 2232e74eeadSLisandro Dalcin #undef __FUNCT__ 2242e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 2252e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 2262e74eeadSLisandro Dalcin { 2272e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 2282e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 2292e74eeadSLisandro Dalcin PetscErrorCode ierr; 2302e74eeadSLisandro Dalcin 2312e74eeadSLisandro Dalcin PetscFunctionBegin; 2322e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 233f23aa3ddSBarry Smith if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 2342e74eeadSLisandro Dalcin #endif 2352e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 2362e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 2372e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 2382e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2392e74eeadSLisandro Dalcin } 2402e74eeadSLisandro Dalcin 2412e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 2422e74eeadSLisandro Dalcin #undef ISG2LMapApply 243b4319ba4SBarry Smith 244b4319ba4SBarry Smith #undef __FUNCT__ 245b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 24613f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 247b4319ba4SBarry Smith { 248dfbe8321SBarry Smith PetscErrorCode ierr; 249b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 250b4319ba4SBarry Smith 251b4319ba4SBarry Smith PetscFunctionBegin; 252b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 253b4319ba4SBarry Smith PetscFunctionReturn(0); 254b4319ba4SBarry Smith } 255b4319ba4SBarry Smith 256b4319ba4SBarry Smith #undef __FUNCT__ 257f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 258f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 259f0006bf2SLisandro Dalcin { 260f0006bf2SLisandro Dalcin PetscErrorCode ierr; 261f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 262f0006bf2SLisandro Dalcin 263f0006bf2SLisandro Dalcin PetscFunctionBegin; 264f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 265f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 266f0006bf2SLisandro Dalcin } 267f0006bf2SLisandro Dalcin 268f0006bf2SLisandro Dalcin #undef __FUNCT__ 2692e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 2702b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2712e74eeadSLisandro Dalcin { 2722e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 2730298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 2742e74eeadSLisandro Dalcin PetscErrorCode ierr; 2752e74eeadSLisandro Dalcin 2762e74eeadSLisandro Dalcin PetscFunctionBegin; 277ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 2782e74eeadSLisandro Dalcin if (n) { 2792e74eeadSLisandro Dalcin ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr); 2802e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 2812e74eeadSLisandro Dalcin } 2822b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 283c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 2842e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2852e74eeadSLisandro Dalcin } 2862e74eeadSLisandro Dalcin 2872e74eeadSLisandro Dalcin #undef __FUNCT__ 288b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 2892b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 290b4319ba4SBarry Smith { 291b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 292dfbe8321SBarry Smith PetscErrorCode ierr; 293f4df32b1SMatthew Knepley PetscInt i; 294b4319ba4SBarry Smith PetscScalar *array; 295b4319ba4SBarry Smith 296b4319ba4SBarry Smith PetscFunctionBegin; 297ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 298b4319ba4SBarry Smith { 299b4319ba4SBarry Smith /* 300b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 301b4319ba4SBarry Smith work properly in the interface nodes. 302b4319ba4SBarry Smith */ 303b4319ba4SBarry Smith Vec counter; 304b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 3050298fd71SBarry Smith ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr); 3062dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 3072dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 308ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 309ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 310ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 311ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3126bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 313b4319ba4SBarry Smith } 314958c9bccSBarry Smith if (!n) { 315b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 316b4319ba4SBarry Smith } else { 317b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 3182205254eSKarl Rupp 319b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 3202b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 321b4319ba4SBarry Smith for (i=0; i<n; i++) { 322f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 323b4319ba4SBarry Smith } 324b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 325b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 326b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 327b4319ba4SBarry Smith } 328b4319ba4SBarry Smith PetscFunctionReturn(0); 329b4319ba4SBarry Smith } 330b4319ba4SBarry Smith 331b4319ba4SBarry Smith #undef __FUNCT__ 332b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 333dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 334b4319ba4SBarry Smith { 335b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 336dfbe8321SBarry Smith PetscErrorCode ierr; 337b4319ba4SBarry Smith 338b4319ba4SBarry Smith PetscFunctionBegin; 339b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 340b4319ba4SBarry Smith PetscFunctionReturn(0); 341b4319ba4SBarry Smith } 342b4319ba4SBarry Smith 343b4319ba4SBarry Smith #undef __FUNCT__ 344b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 345dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 346b4319ba4SBarry Smith { 347b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 348dfbe8321SBarry Smith PetscErrorCode ierr; 349b4319ba4SBarry Smith 350b4319ba4SBarry Smith PetscFunctionBegin; 351b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 352b4319ba4SBarry Smith PetscFunctionReturn(0); 353b4319ba4SBarry Smith } 354b4319ba4SBarry Smith 355b4319ba4SBarry Smith #undef __FUNCT__ 356b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 3577087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 358b4319ba4SBarry Smith { 359b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 360b4319ba4SBarry Smith 361b4319ba4SBarry Smith PetscFunctionBegin; 362b4319ba4SBarry Smith *local = is->A; 363b4319ba4SBarry Smith PetscFunctionReturn(0); 364b4319ba4SBarry Smith } 365b4319ba4SBarry Smith 366b4319ba4SBarry Smith #undef __FUNCT__ 367b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 368b4319ba4SBarry Smith /*@ 369b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 370b4319ba4SBarry Smith 371b4319ba4SBarry Smith Input Parameter: 372b4319ba4SBarry Smith . mat - the matrix 373b4319ba4SBarry Smith 374b4319ba4SBarry Smith Output Parameter: 375b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 376b4319ba4SBarry Smith 377b4319ba4SBarry Smith Level: advanced 378b4319ba4SBarry Smith 379b4319ba4SBarry Smith Notes: 380b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 381b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 382b4319ba4SBarry Smith of the MatSetValues() operation. 383b4319ba4SBarry Smith 384b4319ba4SBarry Smith .seealso: MATIS 385b4319ba4SBarry Smith @*/ 3867087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 387b4319ba4SBarry Smith { 3884ac538c5SBarry Smith PetscErrorCode ierr; 389b4319ba4SBarry Smith 390b4319ba4SBarry Smith PetscFunctionBegin; 3910700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 392b4319ba4SBarry Smith PetscValidPointer(local,2); 3934ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 394b4319ba4SBarry Smith PetscFunctionReturn(0); 395b4319ba4SBarry Smith } 396b4319ba4SBarry Smith 3973b03a366Sstefano_zampini #undef __FUNCT__ 3983b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 3993b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 4003b03a366Sstefano_zampini { 4013b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 4023b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 4033b03a366Sstefano_zampini PetscErrorCode ierr; 4043b03a366Sstefano_zampini 4053b03a366Sstefano_zampini PetscFunctionBegin; 4064e4c7dbeSStefano Zampini if (is->A) { 4073b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 4083b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 409f23aa3ddSBarry Smith if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols); 4104e4c7dbeSStefano Zampini } 4113b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 4123b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 4133b03a366Sstefano_zampini is->A = local; 4143b03a366Sstefano_zampini PetscFunctionReturn(0); 4153b03a366Sstefano_zampini } 4163b03a366Sstefano_zampini 4173b03a366Sstefano_zampini #undef __FUNCT__ 4183b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 4193b03a366Sstefano_zampini /*@ 4203b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 4213b03a366Sstefano_zampini 4223b03a366Sstefano_zampini Input Parameter: 4233b03a366Sstefano_zampini . mat - the matrix 4243b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 4253b03a366Sstefano_zampini 4263b03a366Sstefano_zampini Output Parameter: 4273b03a366Sstefano_zampini 4283b03a366Sstefano_zampini Level: advanced 4293b03a366Sstefano_zampini 4303b03a366Sstefano_zampini Notes: 4313b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 4323b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 4333b03a366Sstefano_zampini 4343b03a366Sstefano_zampini .seealso: MATIS 4353b03a366Sstefano_zampini @*/ 4363b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 4373b03a366Sstefano_zampini { 4383b03a366Sstefano_zampini PetscErrorCode ierr; 4393b03a366Sstefano_zampini 4403b03a366Sstefano_zampini PetscFunctionBegin; 4413b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4423b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 4433b03a366Sstefano_zampini PetscFunctionReturn(0); 4443b03a366Sstefano_zampini } 4453b03a366Sstefano_zampini 4466726f965SBarry Smith #undef __FUNCT__ 4476726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 4486726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 4496726f965SBarry Smith { 4506726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4516726f965SBarry Smith PetscErrorCode ierr; 4526726f965SBarry Smith 4536726f965SBarry Smith PetscFunctionBegin; 4546726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 4556726f965SBarry Smith PetscFunctionReturn(0); 4566726f965SBarry Smith } 4576726f965SBarry Smith 4586726f965SBarry Smith #undef __FUNCT__ 4592e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 4602e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 4612e74eeadSLisandro Dalcin { 4622e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4632e74eeadSLisandro Dalcin PetscErrorCode ierr; 4642e74eeadSLisandro Dalcin 4652e74eeadSLisandro Dalcin PetscFunctionBegin; 4662e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 4672e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4682e74eeadSLisandro Dalcin } 4692e74eeadSLisandro Dalcin 4702e74eeadSLisandro Dalcin #undef __FUNCT__ 4712e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 4722e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 4732e74eeadSLisandro Dalcin { 4742e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4752e74eeadSLisandro Dalcin PetscErrorCode ierr; 4762e74eeadSLisandro Dalcin 4772e74eeadSLisandro Dalcin PetscFunctionBegin; 4782e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 4792e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 4802e74eeadSLisandro Dalcin 4812e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 4822e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 483ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 484ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4852e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4862e74eeadSLisandro Dalcin } 4872e74eeadSLisandro Dalcin 4882e74eeadSLisandro Dalcin #undef __FUNCT__ 4896726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 490ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 4916726f965SBarry Smith { 4926726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4936726f965SBarry Smith PetscErrorCode ierr; 4946726f965SBarry Smith 4956726f965SBarry Smith PetscFunctionBegin; 4964e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 4976726f965SBarry Smith PetscFunctionReturn(0); 4986726f965SBarry Smith } 4996726f965SBarry Smith 500284134d9SBarry Smith #undef __FUNCT__ 501284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 502284134d9SBarry Smith /*@ 503284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 504284134d9SBarry Smith process but not across processes. 505284134d9SBarry Smith 506284134d9SBarry Smith Input Parameters: 507284134d9SBarry Smith + comm - MPI communicator that will share the matrix 5084e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 509284134d9SBarry Smith . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 510284134d9SBarry Smith - map - mapping that defines the global number for each local number 511284134d9SBarry Smith 512284134d9SBarry Smith Output Parameter: 513284134d9SBarry Smith . A - the resulting matrix 514284134d9SBarry Smith 5158e6c10adSSatish Balay Level: advanced 5168e6c10adSSatish Balay 517284134d9SBarry Smith Notes: See MATIS for more details 5188cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 5198cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 5208cda6cd7SBarry Smith plus the ghost points to global indices. 521284134d9SBarry Smith 522284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 523284134d9SBarry Smith @*/ 5244e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 525284134d9SBarry Smith { 526284134d9SBarry Smith PetscErrorCode ierr; 527284134d9SBarry Smith 528284134d9SBarry Smith PetscFunctionBegin; 529284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 530d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 531284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 532284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 5333b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 534784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 535284134d9SBarry Smith PetscFunctionReturn(0); 536284134d9SBarry Smith } 537284134d9SBarry Smith 538b4319ba4SBarry Smith /*MC 539b4319ba4SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 540b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 541b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 542b4319ba4SBarry Smith product is handled "implicitly". 543b4319ba4SBarry Smith 544b4319ba4SBarry Smith Operations Provided: 5456726f965SBarry Smith + MatMult() 5462e74eeadSLisandro Dalcin . MatMultAdd() 5472e74eeadSLisandro Dalcin . MatMultTranspose() 5482e74eeadSLisandro Dalcin . MatMultTransposeAdd() 5496726f965SBarry Smith . MatZeroEntries() 5506726f965SBarry Smith . MatSetOption() 5512e74eeadSLisandro Dalcin . MatZeroRows() 5526726f965SBarry Smith . MatZeroRowsLocal() 5532e74eeadSLisandro Dalcin . MatSetValues() 5546726f965SBarry Smith . MatSetValuesLocal() 5552e74eeadSLisandro Dalcin . MatScale() 5562e74eeadSLisandro Dalcin . MatGetDiagonal() 5576726f965SBarry Smith - MatSetLocalToGlobalMapping() 558b4319ba4SBarry Smith 559b4319ba4SBarry Smith Options Database Keys: 560b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 561b4319ba4SBarry Smith 562b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 563b4319ba4SBarry Smith 564b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 565b4319ba4SBarry Smith 566b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 567b4319ba4SBarry Smith MatISGetLocalMat() 568b4319ba4SBarry Smith 569b4319ba4SBarry Smith Level: advanced 570b4319ba4SBarry Smith 571b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 572b4319ba4SBarry Smith 573b4319ba4SBarry Smith M*/ 574b4319ba4SBarry Smith 575b4319ba4SBarry Smith #undef __FUNCT__ 576b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 5778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 578b4319ba4SBarry Smith { 579dfbe8321SBarry Smith PetscErrorCode ierr; 580b4319ba4SBarry Smith Mat_IS *b; 581b4319ba4SBarry Smith 582b4319ba4SBarry Smith PetscFunctionBegin; 58338f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 584b4319ba4SBarry Smith A->data = (void*)b; 585b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 586b4319ba4SBarry Smith 587b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 5882e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 5892e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 5902e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 591b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 592b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 5932e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 594b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 595f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 5962e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 597b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 598b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 599b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 600b4319ba4SBarry Smith A->ops->view = MatView_IS; 6016726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 6022e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 6032e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 6046726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 60569796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 60669796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 607b4319ba4SBarry Smith 60826283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 60926283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 610b4319ba4SBarry Smith 611b4319ba4SBarry Smith b->A = 0; 612b4319ba4SBarry Smith b->ctx = 0; 613b4319ba4SBarry Smith b->x = 0; 614b4319ba4SBarry Smith b->y = 0; 615bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 616bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 61717667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 618b4319ba4SBarry Smith PetscFunctionReturn(0); 619b4319ba4SBarry Smith } 620