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__ 18ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 19ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 20ad6194a2SStefano Zampini { 21ad6194a2SStefano Zampini PetscErrorCode ierr; 22ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 23ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 24ad6194a2SStefano Zampini Mat B,localmat; 25ad6194a2SStefano Zampini 26ad6194a2SStefano Zampini PetscFunctionBegin; 27ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 28ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 29ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 30ad6194a2SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr); 31ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 32ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 33ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 34ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 35ad6194a2SStefano Zampini *newmat = B; 36ad6194a2SStefano Zampini PetscFunctionReturn(0); 37ad6194a2SStefano Zampini } 38ad6194a2SStefano Zampini 39ad6194a2SStefano Zampini #undef __FUNCT__ 4069796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 4169796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 4269796d55SStefano Zampini { 4369796d55SStefano Zampini PetscErrorCode ierr; 4469796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 4569796d55SStefano Zampini PetscBool local_sym; 4669796d55SStefano Zampini 4769796d55SStefano Zampini PetscFunctionBegin; 4869796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 4969796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 5069796d55SStefano Zampini PetscFunctionReturn(0); 5169796d55SStefano Zampini } 5269796d55SStefano Zampini 5369796d55SStefano Zampini #undef __FUNCT__ 5469796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 5569796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 5669796d55SStefano Zampini { 5769796d55SStefano Zampini PetscErrorCode ierr; 5869796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 5969796d55SStefano Zampini PetscBool local_sym; 6069796d55SStefano Zampini 6169796d55SStefano Zampini PetscFunctionBegin; 6269796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 6369796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 6469796d55SStefano Zampini PetscFunctionReturn(0); 6569796d55SStefano Zampini } 6669796d55SStefano Zampini 6769796d55SStefano Zampini #undef __FUNCT__ 68b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 69dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 70b4319ba4SBarry Smith { 71dfbe8321SBarry Smith PetscErrorCode ierr; 72b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 73b4319ba4SBarry Smith 74b4319ba4SBarry Smith PetscFunctionBegin; 756bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 766bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 776bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 786bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 796bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 80bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 81dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 82bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 83b4319ba4SBarry Smith PetscFunctionReturn(0); 84b4319ba4SBarry Smith } 85b4319ba4SBarry Smith 86b4319ba4SBarry Smith #undef __FUNCT__ 87b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 88dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 89b4319ba4SBarry Smith { 90dfbe8321SBarry Smith PetscErrorCode ierr; 91b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 92b4319ba4SBarry Smith PetscScalar zero = 0.0; 93b4319ba4SBarry Smith 94b4319ba4SBarry Smith PetscFunctionBegin; 95b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 96ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 97ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 98b4319ba4SBarry Smith 99b4319ba4SBarry Smith /* multiply the local matrix */ 100b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 101b4319ba4SBarry Smith 102b4319ba4SBarry Smith /* scatter product back into global memory */ 1032dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 104ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 105ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 106b4319ba4SBarry Smith PetscFunctionReturn(0); 107b4319ba4SBarry Smith } 108b4319ba4SBarry Smith 109b4319ba4SBarry Smith #undef __FUNCT__ 1102e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 1112e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 1122e74eeadSLisandro Dalcin { 113650997f4SStefano Zampini Vec temp_vec; 1142e74eeadSLisandro Dalcin PetscErrorCode ierr; 1152e74eeadSLisandro Dalcin 1162e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 117650997f4SStefano Zampini if (v3 != v2) { 118650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 119650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 120650997f4SStefano Zampini } else { 121650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 122650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 123650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 124650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 125650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 126650997f4SStefano Zampini } 1272e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1282e74eeadSLisandro Dalcin } 1292e74eeadSLisandro Dalcin 1302e74eeadSLisandro Dalcin #undef __FUNCT__ 1312e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 1322e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 1332e74eeadSLisandro Dalcin { 1342e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1352e74eeadSLisandro Dalcin PetscErrorCode ierr; 1362e74eeadSLisandro Dalcin 1372e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 1382e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 139ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 140ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1412e74eeadSLisandro Dalcin 1422e74eeadSLisandro Dalcin /* multiply the local matrix */ 1432e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 1442e74eeadSLisandro Dalcin 1452e74eeadSLisandro Dalcin /* scatter product back into global vector */ 1462e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 147ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 148ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1492e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1502e74eeadSLisandro Dalcin } 1512e74eeadSLisandro Dalcin 1522e74eeadSLisandro Dalcin #undef __FUNCT__ 1532e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1542e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 1552e74eeadSLisandro Dalcin { 156650997f4SStefano Zampini Vec temp_vec; 1572e74eeadSLisandro Dalcin PetscErrorCode ierr; 1582e74eeadSLisandro Dalcin 1592e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 160650997f4SStefano Zampini if (v3 != v2) { 161650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 162650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 163650997f4SStefano Zampini } else { 164650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 165650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 166650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 167650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 168650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 169650997f4SStefano Zampini } 1702e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1712e74eeadSLisandro Dalcin } 1722e74eeadSLisandro Dalcin 1732e74eeadSLisandro Dalcin #undef __FUNCT__ 174b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 175dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 176b4319ba4SBarry Smith { 177b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 178dfbe8321SBarry Smith PetscErrorCode ierr; 179b4319ba4SBarry Smith PetscViewer sviewer; 180b4319ba4SBarry Smith 181b4319ba4SBarry Smith PetscFunctionBegin; 182dd2fa690SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 183b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 184b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 185b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 186b4319ba4SBarry Smith PetscFunctionReturn(0); 187b4319ba4SBarry Smith } 188b4319ba4SBarry Smith 189b4319ba4SBarry Smith #undef __FUNCT__ 190b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 191784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 192b4319ba4SBarry Smith { 193dfbe8321SBarry Smith PetscErrorCode ierr; 1944e4c7dbeSStefano Zampini PetscInt n,bs; 195b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 196b4319ba4SBarry Smith IS from,to; 197b4319ba4SBarry Smith Vec global; 198b4319ba4SBarry Smith 199b4319ba4SBarry Smith PetscFunctionBegin; 200784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 201ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 20270cf5478SStefano Zampini if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 20370cf5478SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 20470cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 20570cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 20670cf5478SStefano Zampini ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 2071c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 20870cf5478SStefano Zampini } 209784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 2106bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 211784ac674SJed Brown is->mapping = rmapping; 212fa7f1dd8SStefano Zampini /* 213fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 214fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 215fa7f1dd8SStefano Zampini */ 216b4319ba4SBarry Smith 217b4319ba4SBarry Smith /* Create the local matrix A */ 218784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 2194e4c7dbeSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 220f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 221f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 2224e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 223ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 224ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 225b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 226b4319ba4SBarry Smith 227b4319ba4SBarry Smith /* Create the local work vectors */ 2284e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 2294e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 2304e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 231ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 232ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 2334e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 234b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 235b4319ba4SBarry Smith 236b4319ba4SBarry Smith /* setup the global to local scatter */ 237b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 238784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 2390298fd71SBarry Smith ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr); 240b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 2416bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 2426bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 2436bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 244b4319ba4SBarry Smith PetscFunctionReturn(0); 245b4319ba4SBarry Smith } 246b4319ba4SBarry Smith 2472e74eeadSLisandro Dalcin #undef __FUNCT__ 2482e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 2492e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 2502e74eeadSLisandro Dalcin { 2512e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 2522e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 2532e74eeadSLisandro Dalcin PetscErrorCode ierr; 2542e74eeadSLisandro Dalcin 2552e74eeadSLisandro Dalcin PetscFunctionBegin; 2562e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 257f23aa3ddSBarry 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); 2582e74eeadSLisandro Dalcin #endif 2592e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 2602e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 2612e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 2622e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2632e74eeadSLisandro Dalcin } 2642e74eeadSLisandro Dalcin 2652e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 2662e74eeadSLisandro Dalcin #undef ISG2LMapApply 267b4319ba4SBarry Smith 268b4319ba4SBarry Smith #undef __FUNCT__ 269b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 27013f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 271b4319ba4SBarry Smith { 272dfbe8321SBarry Smith PetscErrorCode ierr; 273b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 274b4319ba4SBarry Smith 275b4319ba4SBarry Smith PetscFunctionBegin; 276b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 277b4319ba4SBarry Smith PetscFunctionReturn(0); 278b4319ba4SBarry Smith } 279b4319ba4SBarry Smith 280b4319ba4SBarry Smith #undef __FUNCT__ 281f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 282f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 283f0006bf2SLisandro Dalcin { 284f0006bf2SLisandro Dalcin PetscErrorCode ierr; 285f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 286f0006bf2SLisandro Dalcin 287f0006bf2SLisandro Dalcin PetscFunctionBegin; 288f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 289f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 290f0006bf2SLisandro Dalcin } 291f0006bf2SLisandro Dalcin 292f0006bf2SLisandro Dalcin #undef __FUNCT__ 2932e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 2942b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2952e74eeadSLisandro Dalcin { 2962e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 2970298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 2982e74eeadSLisandro Dalcin PetscErrorCode ierr; 2992e74eeadSLisandro Dalcin 3002e74eeadSLisandro Dalcin PetscFunctionBegin; 301ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 3022e74eeadSLisandro Dalcin if (n) { 303785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 3042e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 3052e74eeadSLisandro Dalcin } 3062b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 307c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 3082e74eeadSLisandro Dalcin PetscFunctionReturn(0); 3092e74eeadSLisandro Dalcin } 3102e74eeadSLisandro Dalcin 3112e74eeadSLisandro Dalcin #undef __FUNCT__ 312b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 3132b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 314b4319ba4SBarry Smith { 315b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 316dfbe8321SBarry Smith PetscErrorCode ierr; 317f4df32b1SMatthew Knepley PetscInt i; 318b4319ba4SBarry Smith PetscScalar *array; 319b4319ba4SBarry Smith 320b4319ba4SBarry Smith PetscFunctionBegin; 321ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 322b4319ba4SBarry Smith { 323b4319ba4SBarry Smith /* 324b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 325b4319ba4SBarry Smith work properly in the interface nodes. 326b4319ba4SBarry Smith */ 327b4319ba4SBarry Smith Vec counter; 328b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 3290298fd71SBarry Smith ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr); 3302dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 3312dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 332ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 333ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 334ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 335ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3366bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 337b4319ba4SBarry Smith } 338958c9bccSBarry Smith if (!n) { 339b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 340b4319ba4SBarry Smith } else { 341b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 3422205254eSKarl Rupp 343b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 3442b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 345b4319ba4SBarry Smith for (i=0; i<n; i++) { 346f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 347b4319ba4SBarry Smith } 348b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 349b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 350b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 351b4319ba4SBarry Smith } 352b4319ba4SBarry Smith PetscFunctionReturn(0); 353b4319ba4SBarry Smith } 354b4319ba4SBarry Smith 355b4319ba4SBarry Smith #undef __FUNCT__ 356b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 357dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 358b4319ba4SBarry Smith { 359b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 360dfbe8321SBarry Smith PetscErrorCode ierr; 361b4319ba4SBarry Smith 362b4319ba4SBarry Smith PetscFunctionBegin; 363b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 364b4319ba4SBarry Smith PetscFunctionReturn(0); 365b4319ba4SBarry Smith } 366b4319ba4SBarry Smith 367b4319ba4SBarry Smith #undef __FUNCT__ 368b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 369dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 370b4319ba4SBarry Smith { 371b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 372dfbe8321SBarry Smith PetscErrorCode ierr; 373b4319ba4SBarry Smith 374b4319ba4SBarry Smith PetscFunctionBegin; 375b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 376b4319ba4SBarry Smith PetscFunctionReturn(0); 377b4319ba4SBarry Smith } 378b4319ba4SBarry Smith 379b4319ba4SBarry Smith #undef __FUNCT__ 380b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 3817087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 382b4319ba4SBarry Smith { 383b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 384b4319ba4SBarry Smith 385b4319ba4SBarry Smith PetscFunctionBegin; 386b4319ba4SBarry Smith *local = is->A; 387b4319ba4SBarry Smith PetscFunctionReturn(0); 388b4319ba4SBarry Smith } 389b4319ba4SBarry Smith 390b4319ba4SBarry Smith #undef __FUNCT__ 391b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 392b4319ba4SBarry Smith /*@ 393b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 394b4319ba4SBarry Smith 395b4319ba4SBarry Smith Input Parameter: 396b4319ba4SBarry Smith . mat - the matrix 397b4319ba4SBarry Smith 398b4319ba4SBarry Smith Output Parameter: 399b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 400b4319ba4SBarry Smith 401b4319ba4SBarry Smith Level: advanced 402b4319ba4SBarry Smith 403b4319ba4SBarry Smith Notes: 404b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 405b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 406b4319ba4SBarry Smith of the MatSetValues() operation. 407b4319ba4SBarry Smith 408b4319ba4SBarry Smith .seealso: MATIS 409b4319ba4SBarry Smith @*/ 4107087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 411b4319ba4SBarry Smith { 4124ac538c5SBarry Smith PetscErrorCode ierr; 413b4319ba4SBarry Smith 414b4319ba4SBarry Smith PetscFunctionBegin; 4150700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 416b4319ba4SBarry Smith PetscValidPointer(local,2); 4174ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 418b4319ba4SBarry Smith PetscFunctionReturn(0); 419b4319ba4SBarry Smith } 420b4319ba4SBarry Smith 4213b03a366Sstefano_zampini #undef __FUNCT__ 4223b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 4233b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 4243b03a366Sstefano_zampini { 4253b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 4263b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 4273b03a366Sstefano_zampini PetscErrorCode ierr; 4283b03a366Sstefano_zampini 4293b03a366Sstefano_zampini PetscFunctionBegin; 4304e4c7dbeSStefano Zampini if (is->A) { 4313b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 4323b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 433f23aa3ddSBarry 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); 4344e4c7dbeSStefano Zampini } 4353b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 4363b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 4373b03a366Sstefano_zampini is->A = local; 4383b03a366Sstefano_zampini PetscFunctionReturn(0); 4393b03a366Sstefano_zampini } 4403b03a366Sstefano_zampini 4413b03a366Sstefano_zampini #undef __FUNCT__ 4423b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 4433b03a366Sstefano_zampini /*@ 4443b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 4453b03a366Sstefano_zampini 4463b03a366Sstefano_zampini Input Parameter: 4473b03a366Sstefano_zampini . mat - the matrix 4483b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 4493b03a366Sstefano_zampini 4503b03a366Sstefano_zampini Output Parameter: 4513b03a366Sstefano_zampini 4523b03a366Sstefano_zampini Level: advanced 4533b03a366Sstefano_zampini 4543b03a366Sstefano_zampini Notes: 4553b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 4563b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 4573b03a366Sstefano_zampini 4583b03a366Sstefano_zampini .seealso: MATIS 4593b03a366Sstefano_zampini @*/ 4603b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 4613b03a366Sstefano_zampini { 4623b03a366Sstefano_zampini PetscErrorCode ierr; 4633b03a366Sstefano_zampini 4643b03a366Sstefano_zampini PetscFunctionBegin; 4653b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4663b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 4673b03a366Sstefano_zampini PetscFunctionReturn(0); 4683b03a366Sstefano_zampini } 4693b03a366Sstefano_zampini 4706726f965SBarry Smith #undef __FUNCT__ 4716726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 4726726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 4736726f965SBarry Smith { 4746726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4756726f965SBarry Smith PetscErrorCode ierr; 4766726f965SBarry Smith 4776726f965SBarry Smith PetscFunctionBegin; 4786726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 4796726f965SBarry Smith PetscFunctionReturn(0); 4806726f965SBarry Smith } 4816726f965SBarry Smith 4826726f965SBarry Smith #undef __FUNCT__ 4832e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 4842e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 4852e74eeadSLisandro Dalcin { 4862e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4872e74eeadSLisandro Dalcin PetscErrorCode ierr; 4882e74eeadSLisandro Dalcin 4892e74eeadSLisandro Dalcin PetscFunctionBegin; 4902e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 4912e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4922e74eeadSLisandro Dalcin } 4932e74eeadSLisandro Dalcin 4942e74eeadSLisandro Dalcin #undef __FUNCT__ 4952e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 4962e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 4972e74eeadSLisandro Dalcin { 4982e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4992e74eeadSLisandro Dalcin PetscErrorCode ierr; 5002e74eeadSLisandro Dalcin 5012e74eeadSLisandro Dalcin PetscFunctionBegin; 5022e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 5032e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 5042e74eeadSLisandro Dalcin 5052e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 5062e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 507ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 508ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5092e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5102e74eeadSLisandro Dalcin } 5112e74eeadSLisandro Dalcin 5122e74eeadSLisandro Dalcin #undef __FUNCT__ 5136726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 514ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 5156726f965SBarry Smith { 5166726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 5176726f965SBarry Smith PetscErrorCode ierr; 5186726f965SBarry Smith 5196726f965SBarry Smith PetscFunctionBegin; 5204e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 5216726f965SBarry Smith PetscFunctionReturn(0); 5226726f965SBarry Smith } 5236726f965SBarry Smith 524284134d9SBarry Smith #undef __FUNCT__ 525284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 526284134d9SBarry Smith /*@ 527284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 528284134d9SBarry Smith process but not across processes. 529284134d9SBarry Smith 530284134d9SBarry Smith Input Parameters: 531284134d9SBarry Smith + comm - MPI communicator that will share the matrix 5324e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 533284134d9SBarry Smith . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 534284134d9SBarry Smith - map - mapping that defines the global number for each local number 535284134d9SBarry Smith 536284134d9SBarry Smith Output Parameter: 537284134d9SBarry Smith . A - the resulting matrix 538284134d9SBarry Smith 5398e6c10adSSatish Balay Level: advanced 5408e6c10adSSatish Balay 541284134d9SBarry Smith Notes: See MATIS for more details 5428cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 5438cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 5448cda6cd7SBarry Smith plus the ghost points to global indices. 545284134d9SBarry Smith 546284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 547284134d9SBarry Smith @*/ 5484e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 549284134d9SBarry Smith { 550284134d9SBarry Smith PetscErrorCode ierr; 551284134d9SBarry Smith 552284134d9SBarry Smith PetscFunctionBegin; 553284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 554d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 555284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 556284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 5573b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 558784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 559284134d9SBarry Smith PetscFunctionReturn(0); 560284134d9SBarry Smith } 561284134d9SBarry Smith 562b4319ba4SBarry Smith /*MC 563b4319ba4SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 564b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 565b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 566b4319ba4SBarry Smith product is handled "implicitly". 567b4319ba4SBarry Smith 568b4319ba4SBarry Smith Operations Provided: 5696726f965SBarry Smith + MatMult() 5702e74eeadSLisandro Dalcin . MatMultAdd() 5712e74eeadSLisandro Dalcin . MatMultTranspose() 5722e74eeadSLisandro Dalcin . MatMultTransposeAdd() 5736726f965SBarry Smith . MatZeroEntries() 5746726f965SBarry Smith . MatSetOption() 5752e74eeadSLisandro Dalcin . MatZeroRows() 5766726f965SBarry Smith . MatZeroRowsLocal() 5772e74eeadSLisandro Dalcin . MatSetValues() 5786726f965SBarry Smith . MatSetValuesLocal() 5792e74eeadSLisandro Dalcin . MatScale() 5802e74eeadSLisandro Dalcin . MatGetDiagonal() 5816726f965SBarry Smith - MatSetLocalToGlobalMapping() 582b4319ba4SBarry Smith 583b4319ba4SBarry Smith Options Database Keys: 584b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 585b4319ba4SBarry Smith 586b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 587b4319ba4SBarry Smith 588b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 589b4319ba4SBarry Smith 590b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 591b4319ba4SBarry Smith MatISGetLocalMat() 592b4319ba4SBarry Smith 593b4319ba4SBarry Smith Level: advanced 594b4319ba4SBarry Smith 595b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 596b4319ba4SBarry Smith 597b4319ba4SBarry Smith M*/ 598b4319ba4SBarry Smith 599b4319ba4SBarry Smith #undef __FUNCT__ 600b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 6018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 602b4319ba4SBarry Smith { 603dfbe8321SBarry Smith PetscErrorCode ierr; 604b4319ba4SBarry Smith Mat_IS *b; 605b4319ba4SBarry Smith 606b4319ba4SBarry Smith PetscFunctionBegin; 607b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 608b4319ba4SBarry Smith A->data = (void*)b; 609b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 610b4319ba4SBarry Smith 611b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 6122e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 6132e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 6142e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 615b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 616b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 6172e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 618b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 619f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 6202e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 621b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 622b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 623b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 624b4319ba4SBarry Smith A->ops->view = MatView_IS; 6256726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 6262e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 6272e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 6286726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 62969796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 63069796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 631ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 632b4319ba4SBarry Smith 63326283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 63426283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 635b4319ba4SBarry Smith 636b4319ba4SBarry Smith b->A = 0; 637b4319ba4SBarry Smith b->ctx = 0; 638b4319ba4SBarry Smith b->x = 0; 639b4319ba4SBarry Smith b->y = 0; 640*b09366fdSStefano Zampini b->mapping = 0; 641bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 642bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 64317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 644b4319ba4SBarry Smith PetscFunctionReturn(0); 645b4319ba4SBarry Smith } 646