1b4319ba4SBarry Smith /* 2b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 3b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 4b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 5b4319ba4SBarry Smith product is handled "implicitly". 6b4319ba4SBarry Smith 7b4319ba4SBarry Smith We provide: 8b4319ba4SBarry Smith MatMult() 9b4319ba4SBarry Smith 10b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 11b4319ba4SBarry Smith 12b4319ba4SBarry Smith */ 13b4319ba4SBarry Smith 14b4319ba4SBarry Smith #include "src/mat/impls/is/matis.h" /*I "petscmat.h" I*/ 15b4319ba4SBarry Smith 16b4319ba4SBarry Smith #undef __FUNCT__ 17b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 18dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 19b4319ba4SBarry Smith { 20dfbe8321SBarry Smith PetscErrorCode ierr; 21b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 22b4319ba4SBarry Smith 23b4319ba4SBarry Smith PetscFunctionBegin; 24b4319ba4SBarry Smith if (b->A) { 25b4319ba4SBarry Smith ierr = MatDestroy(b->A);CHKERRQ(ierr); 26b4319ba4SBarry Smith } 27b4319ba4SBarry Smith if (b->ctx) { 28b4319ba4SBarry Smith ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 29b4319ba4SBarry Smith } 30b4319ba4SBarry Smith if (b->x) { 31b4319ba4SBarry Smith ierr = VecDestroy(b->x);CHKERRQ(ierr); 32b4319ba4SBarry Smith } 33b4319ba4SBarry Smith if (b->y) { 34b4319ba4SBarry Smith ierr = VecDestroy(b->y);CHKERRQ(ierr); 35b4319ba4SBarry Smith } 36b4319ba4SBarry Smith if (b->mapping) { 37b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr); 38b4319ba4SBarry Smith } 39b4319ba4SBarry Smith ierr = PetscFree(b);CHKERRQ(ierr); 40901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr); 41b4319ba4SBarry Smith PetscFunctionReturn(0); 42b4319ba4SBarry Smith } 43b4319ba4SBarry Smith 44b4319ba4SBarry Smith #undef __FUNCT__ 45b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 46dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 47b4319ba4SBarry Smith { 48dfbe8321SBarry Smith PetscErrorCode ierr; 49b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 50b4319ba4SBarry Smith PetscScalar zero = 0.0; 51b4319ba4SBarry Smith 52b4319ba4SBarry Smith PetscFunctionBegin; 53b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 54b4319ba4SBarry Smith ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 55b4319ba4SBarry Smith ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 56b4319ba4SBarry Smith 57b4319ba4SBarry Smith /* multiply the local matrix */ 58b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 59b4319ba4SBarry Smith 60b4319ba4SBarry Smith /* scatter product back into global memory */ 61b4319ba4SBarry Smith ierr = VecSet(&zero,y);CHKERRQ(ierr); 62b4319ba4SBarry Smith ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 63b4319ba4SBarry Smith ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 64b4319ba4SBarry Smith 65b4319ba4SBarry Smith PetscFunctionReturn(0); 66b4319ba4SBarry Smith } 67b4319ba4SBarry Smith 68b4319ba4SBarry Smith #undef __FUNCT__ 69b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 70dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 71b4319ba4SBarry Smith { 72b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 73dfbe8321SBarry Smith PetscErrorCode ierr; 74b4319ba4SBarry Smith PetscViewer sviewer; 75b4319ba4SBarry Smith 76b4319ba4SBarry Smith PetscFunctionBegin; 77b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 78b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 79b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 80b4319ba4SBarry Smith PetscFunctionReturn(0); 81b4319ba4SBarry Smith } 82b4319ba4SBarry Smith 83b4319ba4SBarry Smith #undef __FUNCT__ 84b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 85dfbe8321SBarry Smith PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping) 86b4319ba4SBarry Smith { 87dfbe8321SBarry Smith PetscErrorCode ierr; 8813f74950SBarry Smith PetscInt n; 89b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 90b4319ba4SBarry Smith IS from,to; 91b4319ba4SBarry Smith Vec global; 92b4319ba4SBarry Smith 93b4319ba4SBarry Smith PetscFunctionBegin; 94b4319ba4SBarry Smith is->mapping = mapping; 95b4319ba4SBarry Smith ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); 96b4319ba4SBarry Smith 97b4319ba4SBarry Smith /* Create the local matrix A */ 98b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 99b4319ba4SBarry Smith ierr = MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);CHKERRQ(ierr); 100b4319ba4SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr); 101b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 102b4319ba4SBarry Smith 103b4319ba4SBarry Smith /* Create the local work vectors */ 104b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr); 105b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 106b4319ba4SBarry Smith 107b4319ba4SBarry Smith /* setup the global to local scatter */ 108b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 109b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr); 110b4319ba4SBarry Smith ierr = VecCreateMPI(A->comm,A->n,A->N,&global);CHKERRQ(ierr); 111b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 112b4319ba4SBarry Smith ierr = VecDestroy(global);CHKERRQ(ierr); 113b4319ba4SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 114b4319ba4SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 115b4319ba4SBarry Smith PetscFunctionReturn(0); 116b4319ba4SBarry Smith } 117b4319ba4SBarry Smith 118b4319ba4SBarry Smith 119b4319ba4SBarry Smith #undef __FUNCT__ 120b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 12113f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 122b4319ba4SBarry Smith { 123dfbe8321SBarry Smith PetscErrorCode ierr; 124b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 125b4319ba4SBarry Smith 126b4319ba4SBarry Smith PetscFunctionBegin; 127b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 128b4319ba4SBarry Smith PetscFunctionReturn(0); 129b4319ba4SBarry Smith } 130b4319ba4SBarry Smith 131b4319ba4SBarry Smith #undef __FUNCT__ 132b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 133dfbe8321SBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag) 134b4319ba4SBarry Smith { 135b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 136dfbe8321SBarry Smith PetscErrorCode ierr; 13713f74950SBarry Smith PetscInt i,n,*rows; 138b4319ba4SBarry Smith PetscScalar *array; 139b4319ba4SBarry Smith 140b4319ba4SBarry Smith PetscFunctionBegin; 141b4319ba4SBarry Smith { 142b4319ba4SBarry Smith /* 143b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 144b4319ba4SBarry Smith work properly in the interface nodes. 145b4319ba4SBarry Smith */ 146b4319ba4SBarry Smith Vec counter; 147b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 148b4319ba4SBarry Smith ierr = VecCreateMPI(A->comm,A->n,A->N,&counter);CHKERRQ(ierr); 149b4319ba4SBarry Smith ierr = VecSet(&zero,counter);CHKERRQ(ierr); 150b4319ba4SBarry Smith ierr = VecSet(&one,is->x);CHKERRQ(ierr); 151b4319ba4SBarry Smith ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 152b4319ba4SBarry Smith ierr = VecScatterEnd (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 153b4319ba4SBarry Smith ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 154b4319ba4SBarry Smith ierr = VecScatterEnd (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 155b4319ba4SBarry Smith ierr = VecDestroy(counter);CHKERRQ(ierr); 156b4319ba4SBarry Smith } 157b4319ba4SBarry Smith ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr); 158958c9bccSBarry Smith if (!n) { 159b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 160b4319ba4SBarry Smith } else { 161b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 162b4319ba4SBarry Smith ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 163b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 164b4319ba4SBarry Smith ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr); 165b4319ba4SBarry Smith for (i=0; i<n; i++) { 166b4319ba4SBarry Smith ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 167b4319ba4SBarry Smith } 168b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 170b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 171b4319ba4SBarry Smith ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 172b4319ba4SBarry Smith } 173b4319ba4SBarry Smith 174b4319ba4SBarry Smith PetscFunctionReturn(0); 175b4319ba4SBarry Smith } 176b4319ba4SBarry Smith 177b4319ba4SBarry Smith #undef __FUNCT__ 178b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 179dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 180b4319ba4SBarry Smith { 181b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 182dfbe8321SBarry Smith PetscErrorCode ierr; 183b4319ba4SBarry Smith 184b4319ba4SBarry Smith PetscFunctionBegin; 185b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 186b4319ba4SBarry Smith PetscFunctionReturn(0); 187b4319ba4SBarry Smith } 188b4319ba4SBarry Smith 189b4319ba4SBarry Smith #undef __FUNCT__ 190b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 191dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 192b4319ba4SBarry Smith { 193b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 194dfbe8321SBarry Smith PetscErrorCode ierr; 195b4319ba4SBarry Smith 196b4319ba4SBarry Smith PetscFunctionBegin; 197b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 198b4319ba4SBarry Smith PetscFunctionReturn(0); 199b4319ba4SBarry Smith } 200b4319ba4SBarry Smith 201b4319ba4SBarry Smith EXTERN_C_BEGIN 202b4319ba4SBarry Smith #undef __FUNCT__ 203b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 204dfbe8321SBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 205b4319ba4SBarry Smith { 206b4319ba4SBarry Smith Mat_IS *is = (Mat_IS *)mat->data; 207b4319ba4SBarry Smith 208b4319ba4SBarry Smith PetscFunctionBegin; 209b4319ba4SBarry Smith *local = is->A; 210b4319ba4SBarry Smith PetscFunctionReturn(0); 211b4319ba4SBarry Smith } 212b4319ba4SBarry Smith EXTERN_C_END 213b4319ba4SBarry Smith 214b4319ba4SBarry Smith #undef __FUNCT__ 215b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 216b4319ba4SBarry Smith /*@ 217b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 218b4319ba4SBarry Smith 219b4319ba4SBarry Smith Input Parameter: 220b4319ba4SBarry Smith . mat - the matrix 221b4319ba4SBarry Smith 222b4319ba4SBarry Smith Output Parameter: 223b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 224b4319ba4SBarry Smith 225b4319ba4SBarry Smith Level: advanced 226b4319ba4SBarry Smith 227b4319ba4SBarry Smith Notes: 228b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 229b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 230b4319ba4SBarry Smith of the MatSetValues() operation. 231b4319ba4SBarry Smith 232b4319ba4SBarry Smith .seealso: MATIS 233b4319ba4SBarry Smith @*/ 234dfbe8321SBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 235b4319ba4SBarry Smith { 236dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,Mat *); 237b4319ba4SBarry Smith 238b4319ba4SBarry Smith PetscFunctionBegin; 239b4319ba4SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 240b4319ba4SBarry Smith PetscValidPointer(local,2); 241b4319ba4SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr); 242b4319ba4SBarry Smith if (f) { 243b4319ba4SBarry Smith ierr = (*f)(mat,local);CHKERRQ(ierr); 244b4319ba4SBarry Smith } else { 245b4319ba4SBarry Smith local = 0; 246b4319ba4SBarry Smith } 247b4319ba4SBarry Smith PetscFunctionReturn(0); 248b4319ba4SBarry Smith } 249b4319ba4SBarry Smith 250*6726f965SBarry Smith #undef __FUNCT__ 251*6726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 252*6726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 253*6726f965SBarry Smith { 254*6726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 255*6726f965SBarry Smith PetscErrorCode ierr; 256*6726f965SBarry Smith 257*6726f965SBarry Smith PetscFunctionBegin; 258*6726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 259*6726f965SBarry Smith PetscFunctionReturn(0); 260*6726f965SBarry Smith } 261*6726f965SBarry Smith 262*6726f965SBarry Smith #undef __FUNCT__ 263*6726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 264*6726f965SBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op) 265*6726f965SBarry Smith { 266*6726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 267*6726f965SBarry Smith PetscErrorCode ierr; 268*6726f965SBarry Smith 269*6726f965SBarry Smith PetscFunctionBegin; 270*6726f965SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 271*6726f965SBarry Smith PetscFunctionReturn(0); 272*6726f965SBarry Smith } 273*6726f965SBarry Smith 274b4319ba4SBarry Smith /*MC 275b4319ba4SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 276b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 277b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 278b4319ba4SBarry Smith product is handled "implicitly". 279b4319ba4SBarry Smith 280b4319ba4SBarry Smith Operations Provided: 281*6726f965SBarry Smith + MatMult() 282*6726f965SBarry Smith . MatZeroEntries() 283*6726f965SBarry Smith . MatSetOption() 284*6726f965SBarry Smith . MatZeroRowsLocal() 285*6726f965SBarry Smith . MatSetValuesLocal() 286*6726f965SBarry Smith - MatSetLocalToGlobalMapping() 287b4319ba4SBarry Smith 288b4319ba4SBarry Smith Options Database Keys: 289b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 290b4319ba4SBarry Smith 291b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 292b4319ba4SBarry Smith 293b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 294b4319ba4SBarry Smith 295b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 296b4319ba4SBarry Smith MatISGetLocalMat() 297b4319ba4SBarry Smith 298b4319ba4SBarry Smith Level: advanced 299b4319ba4SBarry Smith 300b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 301b4319ba4SBarry Smith 302b4319ba4SBarry Smith M*/ 303b4319ba4SBarry Smith 304b4319ba4SBarry Smith EXTERN_C_BEGIN 305b4319ba4SBarry Smith #undef __FUNCT__ 306b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 307dfbe8321SBarry Smith PetscErrorCode MatCreate_IS(Mat A) 308b4319ba4SBarry Smith { 309dfbe8321SBarry Smith PetscErrorCode ierr; 310b4319ba4SBarry Smith Mat_IS *b; 311b4319ba4SBarry Smith 312b4319ba4SBarry Smith PetscFunctionBegin; 313b4319ba4SBarry Smith ierr = PetscNew(Mat_IS,&b);CHKERRQ(ierr); 314b4319ba4SBarry Smith A->data = (void*)b; 315b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 316b4319ba4SBarry Smith A->factor = 0; 317b4319ba4SBarry Smith A->mapping = 0; 318b4319ba4SBarry Smith 319b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 320b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 321b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 322b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 323b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 324b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 325b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 326b4319ba4SBarry Smith A->ops->view = MatView_IS; 327*6726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 328*6726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 329b4319ba4SBarry Smith 330b4319ba4SBarry Smith ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr); 331b4319ba4SBarry Smith ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr); 33213f74950SBarry Smith ierr = MPI_Scan(&A->m,&b->rend,1,MPIU_INT,MPI_SUM,A->comm);CHKERRQ(ierr); 333b4319ba4SBarry Smith b->rstart = b->rend - A->m; 334b4319ba4SBarry Smith 335b4319ba4SBarry Smith b->A = 0; 336b4319ba4SBarry Smith b->ctx = 0; 337b4319ba4SBarry Smith b->x = 0; 338b4319ba4SBarry Smith b->y = 0; 339b4319ba4SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 340b4319ba4SBarry Smith 341b4319ba4SBarry Smith PetscFunctionReturn(0); 342b4319ba4SBarry Smith } 343b4319ba4SBarry Smith EXTERN_C_END 344b4319ba4SBarry Smith 345b4319ba4SBarry Smith 346b4319ba4SBarry Smith 347b4319ba4SBarry Smith 348b4319ba4SBarry Smith 349b4319ba4SBarry Smith 350b4319ba4SBarry Smith 351b4319ba4SBarry Smith 352b4319ba4SBarry Smith 353b4319ba4SBarry Smith 354b4319ba4SBarry Smith 355b4319ba4SBarry Smith 356b4319ba4SBarry Smith 357