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