#define PETSCMAT_DLL /* Creates a matrix class for using the Neumann-Neumann type preconditioners. This stores the matrices in globally unassembled form. Each processor assembles only its local Neumann problem and the parallel matrix vector product is handled "implicitly". We provide: MatMult() Currently this allows for only one subdomain per processor. */ #include "src/mat/impls/is/matis.h" /*I "petscmat.h" I*/ #undef __FUNCT__ #define __FUNCT__ "MatDestroy_IS" PetscErrorCode MatDestroy_IS(Mat A) { PetscErrorCode ierr; Mat_IS *b = (Mat_IS*)A->data; PetscFunctionBegin; if (b->A) { ierr = MatDestroy(b->A);CHKERRQ(ierr); } if (b->ctx) { ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); } if (b->x) { ierr = VecDestroy(b->x);CHKERRQ(ierr); } if (b->y) { ierr = VecDestroy(b->y);CHKERRQ(ierr); } if (b->mapping) { ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr); } ierr = PetscFree(b);CHKERRQ(ierr); ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_IS" PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) { PetscErrorCode ierr; Mat_IS *is = (Mat_IS*)A->data; PetscScalar zero = 0.0; PetscFunctionBegin; /* scatter the global vector x into the local work vector */ ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); /* multiply the local matrix */ ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); /* scatter product back into global memory */ ierr = VecSet(y,zero);CHKERRQ(ierr); ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatView_IS" PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) { Mat_IS *a = (Mat_IS*)A->data; PetscErrorCode ierr; PetscViewer sviewer; PetscFunctionBegin; ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); ierr = MatView(a->A,sviewer);CHKERRQ(ierr); ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping) { PetscErrorCode ierr; PetscInt n; Mat_IS *is = (Mat_IS*)A->data; IS from,to; Vec global; PetscFunctionBegin; is->mapping = mapping; ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); /* Create the local matrix A */ ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr); ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); /* Create the local work vectors */ ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr); ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); /* setup the global to local scatter */ ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr); ierr = VecCreateMPI(A->comm,A->cmap.n,A->cmap.N,&global);CHKERRQ(ierr); ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); ierr = VecDestroy(global);CHKERRQ(ierr); ierr = ISDestroy(to);CHKERRQ(ierr); ierr = ISDestroy(from);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesLocal_IS" PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) { PetscErrorCode ierr; Mat_IS *is = (Mat_IS*)A->data; PetscFunctionBegin; ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatZeroRowsLocal_IS" PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag) { Mat_IS *is = (Mat_IS*)A->data; PetscErrorCode ierr; PetscInt i; PetscScalar *array; PetscFunctionBegin; { /* Set up is->x as a "counting vector". This is in order to MatMult_IS work properly in the interface nodes. */ Vec counter; PetscScalar one=1.0, zero=0.0; ierr = VecCreateMPI(A->comm,A->cmap.n,A->cmap.N,&counter);CHKERRQ(ierr); ierr = VecSet(counter,zero);CHKERRQ(ierr); ierr = VecSet(is->x,one);CHKERRQ(ierr); ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); ierr = VecScatterEnd (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); ierr = VecScatterEnd (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); ierr = VecDestroy(counter);CHKERRQ(ierr); } if (!n) { is->pure_neumann = PETSC_TRUE; } else { is->pure_neumann = PETSC_FALSE; ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); ierr = MatZeroRows(is->A,n,rows,diag);CHKERRQ(ierr); for (i=0; iA,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); } ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatAssemblyBegin_IS" PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) { Mat_IS *is = (Mat_IS*)A->data; PetscErrorCode ierr; PetscFunctionBegin; ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatAssemblyEnd_IS" PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) { Mat_IS *is = (Mat_IS*)A->data; PetscErrorCode ierr; PetscFunctionBegin; ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatISGetLocalMat_IS" PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat_IS(Mat mat,Mat *local) { Mat_IS *is = (Mat_IS *)mat->data; PetscFunctionBegin; *local = is->A; PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatISGetLocalMat" /*@ MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. Input Parameter: . mat - the matrix Output Parameter: . local - the local matrix usually MATSEQAIJ Level: advanced Notes: This can be called if you have precomputed the nonzero structure of the matrix and want to provide it to the inner matrix object to improve the performance of the MatSetValues() operation. .seealso: MATIS @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat mat,Mat *local) { PetscErrorCode ierr,(*f)(Mat,Mat *); PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidPointer(local,2); ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr); if (f) { ierr = (*f)(mat,local);CHKERRQ(ierr); } else { local = 0; } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatZeroEntries_IS" PetscErrorCode MatZeroEntries_IS(Mat A) { Mat_IS *a = (Mat_IS*)A->data; PetscErrorCode ierr; PetscFunctionBegin; ierr = MatZeroEntries(a->A);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetOption_IS" PetscErrorCode MatSetOption_IS(Mat A,MatOption op) { Mat_IS *a = (Mat_IS*)A->data; PetscErrorCode ierr; PetscFunctionBegin; ierr = MatSetOption(a->A,op);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCreateIS" /*@ MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each process but not across processes. Input Parameters: + comm - MPI communicator that will share the matrix . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products - map - mapping that defines the global number for each local number Output Parameter: . A - the resulting matrix Level: advanced Notes: See MATIS for more details m and n are NOT related to the size of the map .seealso: MATIS, MatSetLocalToGlobalMapping() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) { PetscErrorCode ierr; PetscFunctionBegin; ierr = MatCreate(comm,A);CHKERRQ(ierr); ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); ierr = MatSetLocalToGlobalMapping(*A,map);CHKERRQ(ierr); PetscFunctionReturn(0); } /*MC MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. This stores the matrices in globally unassembled form. Each processor assembles only its local Neumann problem and the parallel matrix vector product is handled "implicitly". Operations Provided: + MatMult() . MatZeroEntries() . MatSetOption() . MatZeroRowsLocal() . MatSetValuesLocal() - MatSetLocalToGlobalMapping() Options Database Keys: . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() Notes: Options prefix for the inner matrix are given by -is_mat_xxx You must call MatSetLocalToGlobalMapping() before using this matrix type. You can do matrix preallocation on the local matrix after you obtain it with MatISGetLocalMat() Level: advanced .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() M*/ EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatCreate_IS" PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_IS(Mat A) { PetscErrorCode ierr; Mat_IS *b; PetscFunctionBegin; ierr = PetscNew(Mat_IS,&b);CHKERRQ(ierr); A->data = (void*)b; ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); A->factor = 0; A->mapping = 0; A->ops->mult = MatMult_IS; A->ops->destroy = MatDestroy_IS; A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; A->ops->setvalueslocal = MatSetValuesLocal_IS; A->ops->zerorowslocal = MatZeroRowsLocal_IS; A->ops->assemblybegin = MatAssemblyBegin_IS; A->ops->assemblyend = MatAssemblyEnd_IS; A->ops->view = MatView_IS; A->ops->zeroentries = MatZeroEntries_IS; A->ops->setoption = MatSetOption_IS; ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr); ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr); b->A = 0; b->ctx = 0; b->x = 0; b->y = 0; ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_END