Actual source code: matis.c
1: /*
2: Creates a matrix class for using the Neumann-Neumann type preconditioners.
3: This stores the matrices in globally unassembled form. Each processor
4: assembles only its local Neumann problem and the parallel matrix vector
5: product is handled "implicitly".
7: We provide:
8: MatMult()
10: Currently this allows for only one subdomain per processor.
12: */
14: #include src/mat/impls/is/matis.h
18: PetscErrorCode MatDestroy_IS(Mat A)
19: {
21: Mat_IS *b = (Mat_IS*)A->data;
24: if (b->A) {
25: MatDestroy(b->A);
26: }
27: if (b->ctx) {
28: VecScatterDestroy(b->ctx);
29: }
30: if (b->x) {
31: VecDestroy(b->x);
32: }
33: if (b->y) {
34: VecDestroy(b->y);
35: }
36: if (b->mapping) {
37: ISLocalToGlobalMappingDestroy(b->mapping);
38: }
39: PetscFree(b);
40: PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);
41: return(0);
42: }
46: PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
47: {
49: Mat_IS *is = (Mat_IS*)A->data;
50: PetscScalar zero = 0.0;
53: /* scatter the global vector x into the local work vector */
54: VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
55: VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
57: /* multiply the local matrix */
58: MatMult(is->A,is->x,is->y);
60: /* scatter product back into global memory */
61: VecSet(&zero,y);
62: VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);
63: VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);
65: return(0);
66: }
70: PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
71: {
72: Mat_IS *a = (Mat_IS*)A->data;
74: PetscViewer sviewer;
77: PetscViewerGetSingleton(viewer,&sviewer);
78: MatView(a->A,sviewer);
79: PetscViewerRestoreSingleton(viewer,&sviewer);
80: return(0);
81: }
85: PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
86: {
88: PetscInt n;
89: Mat_IS *is = (Mat_IS*)A->data;
90: IS from,to;
91: Vec global;
94: is->mapping = mapping;
95: PetscObjectReference((PetscObject)mapping);
97: /* Create the local matrix A */
98: ISLocalToGlobalMappingGetSize(mapping,&n);
99: MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);
100: PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");
101: MatSetFromOptions(is->A);
103: /* Create the local work vectors */
104: VecCreateSeq(PETSC_COMM_SELF,n,&is->x);
105: VecDuplicate(is->x,&is->y);
107: /* setup the global to local scatter */
108: ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);
109: ISLocalToGlobalMappingApplyIS(mapping,to,&from);
110: VecCreateMPI(A->comm,A->n,A->N,&global);
111: VecScatterCreate(global,from,is->x,to,&is->ctx);
112: VecDestroy(global);
113: ISDestroy(to);
114: ISDestroy(from);
115: return(0);
116: }
121: PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
122: {
124: Mat_IS *is = (Mat_IS*)A->data;
127: MatSetValues(is->A,m,rows,n,cols,values,addv);
128: return(0);
129: }
133: PetscErrorCode MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag)
134: {
135: Mat_IS *is = (Mat_IS*)A->data;
137: PetscInt i,n,*rows;
138: PetscScalar *array;
141: {
142: /*
143: Set up is->x as a "counting vector". This is in order to MatMult_IS
144: work properly in the interface nodes.
145: */
146: Vec counter;
147: PetscScalar one=1.0, zero=0.0;
148: VecCreateMPI(A->comm,A->n,A->N,&counter);
149: VecSet(&zero,counter);
150: VecSet(&one,is->x);
151: VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);
152: VecScatterEnd (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);
153: VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
154: VecScatterEnd (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
155: VecDestroy(counter);
156: }
157: ISGetLocalSize(isrows,&n);
158: if (!n) {
159: is->pure_neumann = PETSC_TRUE;
160: } else {
161: is->pure_neumann = PETSC_FALSE;
162: ISGetIndices(isrows,&rows);
163: VecGetArray(is->x,&array);
164: MatZeroRows(is->A,isrows,diag);
165: for (i=0; i<n; i++) {
166: MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);
167: }
168: MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);
169: MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);
170: VecRestoreArray(is->x,&array);
171: ISRestoreIndices(isrows,&rows);
172: }
174: return(0);
175: }
179: PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
180: {
181: Mat_IS *is = (Mat_IS*)A->data;
185: MatAssemblyBegin(is->A,type);
186: return(0);
187: }
191: PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
192: {
193: Mat_IS *is = (Mat_IS*)A->data;
197: MatAssemblyEnd(is->A,type);
198: return(0);
199: }
204: PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
205: {
206: Mat_IS *is = (Mat_IS *)mat->data;
207:
209: *local = is->A;
210: return(0);
211: }
216: /*@
217: MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
219: Input Parameter:
220: . mat - the matrix
222: Output Parameter:
223: . local - the local matrix usually MATSEQAIJ
225: Level: advanced
227: Notes:
228: This can be called if you have precomputed the nonzero structure of the
229: matrix and want to provide it to the inner matrix object to improve the performance
230: of the MatSetValues() operation.
232: .seealso: MATIS
233: @*/
234: PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
235: {
236: PetscErrorCode ierr,(*f)(Mat,Mat *);
241: PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);
242: if (f) {
243: (*f)(mat,local);
244: } else {
245: local = 0;
246: }
247: return(0);
248: }
250: /*MC
251: MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
252: This stores the matrices in globally unassembled form. Each processor
253: assembles only its local Neumann problem and the parallel matrix vector
254: product is handled "implicitly".
256: Operations Provided:
257: . MatMult
259: Options Database Keys:
260: . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
262: Notes: Options prefix for the inner matrix are given by -is_mat_xxx
263:
264: You must call MatSetLocalToGlobalMapping() before using this matrix type.
266: You can do matrix preallocation on the local matrix after you obtain it with
267: MatISGetLocalMat()
269: Level: advanced
271: .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
273: M*/
278: PetscErrorCode MatCreate_IS(Mat A)
279: {
281: Mat_IS *b;
284: PetscNew(Mat_IS,&b);
285: A->data = (void*)b;
286: PetscMemzero(b,sizeof(Mat_IS));
287: PetscMemzero(A->ops,sizeof(struct _MatOps));
288: A->factor = 0;
289: A->mapping = 0;
291: A->ops->mult = MatMult_IS;
292: A->ops->destroy = MatDestroy_IS;
293: A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
294: A->ops->setvalueslocal = MatSetValuesLocal_IS;
295: A->ops->zerorowslocal = MatZeroRowsLocal_IS;
296: A->ops->assemblybegin = MatAssemblyBegin_IS;
297: A->ops->assemblyend = MatAssemblyEnd_IS;
298: A->ops->view = MatView_IS;
300: PetscSplitOwnership(A->comm,&A->m,&A->M);
301: PetscSplitOwnership(A->comm,&A->n,&A->N);
302: MPI_Scan(&A->m,&b->rend,1,MPIU_INT,MPI_SUM,A->comm);
303: b->rstart = b->rend - A->m;
305: b->A = 0;
306: b->ctx = 0;
307: b->x = 0;
308: b->y = 0;
309: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);
311: return(0);
312: }