Actual source code: matis.c
1: #define PETSCMAT_DLL
3: /*
4: Creates a matrix class for using the Neumann-Neumann type preconditioners.
5: This stores the matrices in globally unassembled form. Each processor
6: assembles only its local Neumann problem and the parallel matrix vector
7: product is handled "implicitly".
9: We provide:
10: MatMult()
12: Currently this allows for only one subdomain per processor.
14: */
16: #include src/mat/impls/is/matis.h
20: PetscErrorCode MatDestroy_IS(Mat A)
21: {
23: Mat_IS *b = (Mat_IS*)A->data;
26: if (b->A) {
27: MatDestroy(b->A);
28: }
29: if (b->ctx) {
30: VecScatterDestroy(b->ctx);
31: }
32: if (b->x) {
33: VecDestroy(b->x);
34: }
35: if (b->y) {
36: VecDestroy(b->y);
37: }
38: if (b->mapping) {
39: ISLocalToGlobalMappingDestroy(b->mapping);
40: }
41: PetscFree(b);
42: PetscObjectChangeTypeName((PetscObject)A,0);
43: PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);
44: return(0);
45: }
49: PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
50: {
52: Mat_IS *is = (Mat_IS*)A->data;
53: PetscScalar zero = 0.0;
56: /* scatter the global vector x into the local work vector */
57: VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
58: VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
60: /* multiply the local matrix */
61: MatMult(is->A,is->x,is->y);
63: /* scatter product back into global memory */
64: VecSet(y,zero);
65: VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
66: VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
68: return(0);
69: }
73: static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
74: {
76: Mat_IS *is = (Mat_IS*)A->data;
79: /* scatter the global vector v1 into the local work vector */
80: VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
81: VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
82: VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
83: VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
85: /* multiply the local matrix */
86: MatMultAdd(is->A,is->x,is->y,is->y);
88: /* scatter result back into global vector */
89: VecSet(v3,0);
90: VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
91: VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
92: return(0);
93: }
97: PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
98: {
99: Mat_IS *is = (Mat_IS*)A->data;
103: /* scatter the global vector x into the local work vector */
104: VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
105: VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
107: /* multiply the local matrix */
108: MatMultTranspose(is->A,is->x,is->y);
110: /* scatter product back into global vector */
111: VecSet(y,0);
112: VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
113: VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
114: return(0);
115: }
119: PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
120: {
121: Mat_IS *is = (Mat_IS*)A->data;
125: /* scatter the global vectors v1 and v2 into the local work vectors */
126: VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
127: VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
128: VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
129: VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
131: /* multiply the local matrix */
132: MatMultTransposeAdd(is->A,is->x,is->y,is->y);
134: /* scatter result back into global vector */
135: VecSet(v3,0);
136: VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
137: VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
138: return(0);
139: }
143: PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
144: {
145: Mat_IS *a = (Mat_IS*)A->data;
147: PetscViewer sviewer;
150: PetscViewerGetSingleton(viewer,&sviewer);
151: MatView(a->A,sviewer);
152: PetscViewerRestoreSingleton(viewer,&sviewer);
153: return(0);
154: }
158: PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
159: {
161: PetscInt n;
162: Mat_IS *is = (Mat_IS*)A->data;
163: IS from,to;
164: Vec global;
167: if (is->mapping) {
168: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
169: }
171: PetscObjectReference((PetscObject)mapping);
172: if (is->mapping) { ISLocalToGlobalMappingDestroy(is->mapping); }
173: is->mapping = mapping;
175: /* Create the local matrix A */
176: ISLocalToGlobalMappingGetSize(mapping,&n);
177: MatCreate(PETSC_COMM_SELF,&is->A);
178: MatSetSizes(is->A,n,n,n,n);
179: MatSetOptionsPrefix(is->A,"is");
180: MatSetFromOptions(is->A);
182: /* Create the local work vectors */
183: VecCreateSeq(PETSC_COMM_SELF,n,&is->x);
184: VecDuplicate(is->x,&is->y);
186: /* setup the global to local scatter */
187: ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);
188: ISLocalToGlobalMappingApplyIS(mapping,to,&from);
189: VecCreateMPI(A->comm,A->cmap.n,A->cmap.N,&global);
190: VecScatterCreate(global,from,is->x,to,&is->ctx);
191: VecDestroy(global);
192: ISDestroy(to);
193: ISDestroy(from);
194: return(0);
195: }
197: #define ISG2LMapApply(mapping,n,in,out) 0;\
198: if (!(mapping)->globals) {\
199: PetscErrorCode _ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\
200: }\
201: {\
202: PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\
203: for (_i=0; _i<n; _i++) {\
204: if (in[_i] < 0) out[_i] = in[_i];\
205: else if (in[_i] < _start) out[_i] = -1;\
206: else if (in[_i] > _end) out[_i] = -1;\
207: else out[_i] = _globals[in[_i] - _start];\
208: }\
209: }
213: PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
214: {
215: Mat_IS *is = (Mat_IS*)mat->data;
216: PetscInt rows_l[2048],cols_l[2048];
220: #if defined(PETSC_USE_DEBUG)
221: if (m > 2048 || n > 2048) {
222: SETERRQ2(PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
223: }
224: #endif
225: ISG2LMapApply(is->mapping,m,rows,rows_l);
226: ISG2LMapApply(is->mapping,n,cols,cols_l);
227: MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);
228: return(0);
229: }
231: #undef ISG2LMapSetUp
232: #undef ISG2LMapApply
236: PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
237: {
239: Mat_IS *is = (Mat_IS*)A->data;
242: MatSetValues(is->A,m,rows,n,cols,values,addv);
243: return(0);
244: }
248: PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
249: {
250: Mat_IS *is = (Mat_IS*)A->data;
251: PetscInt n_l=0, *rows_l = PETSC_NULL;
255: if (n) {
256: PetscMalloc(n*sizeof(PetscInt),&rows_l);
257: ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);
258: }
259: MatZeroRowsLocal(A,n_l,rows_l,diag);
260: if (rows_l) { PetscFree(rows_l); }
261: return(0);
262: }
266: PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
267: {
268: Mat_IS *is = (Mat_IS*)A->data;
270: PetscInt i;
271: PetscScalar *array;
274: {
275: /*
276: Set up is->x as a "counting vector". This is in order to MatMult_IS
277: work properly in the interface nodes.
278: */
279: Vec counter;
280: PetscScalar one=1.0, zero=0.0;
281: VecCreateMPI(A->comm,A->cmap.n,A->cmap.N,&counter);
282: VecSet(counter,zero);
283: VecSet(is->x,one);
284: VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);
285: VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);
286: VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);
287: VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);
288: VecDestroy(counter);
289: }
290: if (!n) {
291: is->pure_neumann = PETSC_TRUE;
292: } else {
293: is->pure_neumann = PETSC_FALSE;
294: VecGetArray(is->x,&array);
295: MatZeroRows(is->A,n,rows,diag);
296: for (i=0; i<n; i++) {
297: MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);
298: }
299: MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);
300: MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);
301: VecRestoreArray(is->x,&array);
302: }
304: return(0);
305: }
309: PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
310: {
311: Mat_IS *is = (Mat_IS*)A->data;
315: MatAssemblyBegin(is->A,type);
316: return(0);
317: }
321: PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
322: {
323: Mat_IS *is = (Mat_IS*)A->data;
327: MatAssemblyEnd(is->A,type);
328: return(0);
329: }
334: PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
335: {
336: Mat_IS *is = (Mat_IS *)mat->data;
337:
339: *local = is->A;
340: return(0);
341: }
346: /*@
347: MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
349: Input Parameter:
350: . mat - the matrix
352: Output Parameter:
353: . local - the local matrix usually MATSEQAIJ
355: Level: advanced
357: Notes:
358: This can be called if you have precomputed the nonzero structure of the
359: matrix and want to provide it to the inner matrix object to improve the performance
360: of the MatSetValues() operation.
362: .seealso: MATIS
363: @*/
364: PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
365: {
366: PetscErrorCode ierr,(*f)(Mat,Mat *);
371: PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);
372: if (f) {
373: (*f)(mat,local);
374: } else {
375: local = 0;
376: }
377: return(0);
378: }
382: PetscErrorCode MatZeroEntries_IS(Mat A)
383: {
384: Mat_IS *a = (Mat_IS*)A->data;
388: MatZeroEntries(a->A);
389: return(0);
390: }
394: PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
395: {
396: Mat_IS *is = (Mat_IS*)A->data;
400: MatScale(is->A,a);
401: return(0);
402: }
406: PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
407: {
408: Mat_IS *is = (Mat_IS*)A->data;
412: /* get diagonal of the local matrix */
413: MatGetDiagonal(is->A,is->x);
415: /* scatter diagonal back into global vector */
416: VecSet(v,0);
417: VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);
418: VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);
419: return(0);
420: }
424: PetscErrorCode MatSetOption_IS(Mat A,MatOption op)
425: {
426: Mat_IS *a = (Mat_IS*)A->data;
430: MatSetOption(a->A,op);
431: return(0);
432: }
436: /*@
437: MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
438: process but not across processes.
440: Input Parameters:
441: + comm - MPI communicator that will share the matrix
442: . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
443: - map - mapping that defines the global number for each local number
445: Output Parameter:
446: . A - the resulting matrix
448: Level: advanced
450: Notes: See MATIS for more details
451: m and n are NOT related to the size of the map
453: .seealso: MATIS, MatSetLocalToGlobalMapping()
454: @*/
455: PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
456: {
460: MatCreate(comm,A);
461: MatSetSizes(*A,m,n,M,N);
462: MatSetType(*A,MATIS);
463: MatSetLocalToGlobalMapping(*A,map);
464: return(0);
465: }
467: /*MC
468: MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
469: This stores the matrices in globally unassembled form. Each processor
470: assembles only its local Neumann problem and the parallel matrix vector
471: product is handled "implicitly".
473: Operations Provided:
474: + MatMult()
475: . MatMultAdd()
476: . MatMultTranspose()
477: . MatMultTransposeAdd()
478: . MatZeroEntries()
479: . MatSetOption()
480: . MatZeroRows()
481: . MatZeroRowsLocal()
482: . MatSetValues()
483: . MatSetValuesLocal()
484: . MatScale()
485: . MatGetDiagonal()
486: - MatSetLocalToGlobalMapping()
488: Options Database Keys:
489: . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
491: Notes: Options prefix for the inner matrix are given by -is_mat_xxx
492:
493: You must call MatSetLocalToGlobalMapping() before using this matrix type.
495: You can do matrix preallocation on the local matrix after you obtain it with
496: MatISGetLocalMat()
498: Level: advanced
500: .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
502: M*/
507: PetscErrorCode MatCreate_IS(Mat A)
508: {
510: Mat_IS *b;
513: PetscNew(Mat_IS,&b);
514: A->data = (void*)b;
515: PetscMemzero(A->ops,sizeof(struct _MatOps));
516: A->factor = 0;
517: A->mapping = 0;
519: A->ops->mult = MatMult_IS;
520: A->ops->multadd = MatMultAdd_IS;
521: A->ops->multtranspose = MatMultTranspose_IS;
522: A->ops->multtransposeadd = MatMultTransposeAdd_IS;
523: A->ops->destroy = MatDestroy_IS;
524: A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
525: A->ops->setvalues = MatSetValues_IS;
526: A->ops->setvalueslocal = MatSetValuesLocal_IS;
527: A->ops->zerorows = MatZeroRows_IS;
528: A->ops->zerorowslocal = MatZeroRowsLocal_IS;
529: A->ops->assemblybegin = MatAssemblyBegin_IS;
530: A->ops->assemblyend = MatAssemblyEnd_IS;
531: A->ops->view = MatView_IS;
532: A->ops->zeroentries = MatZeroEntries_IS;
533: A->ops->scale = MatScale_IS;
534: A->ops->getdiagonal = MatGetDiagonal_IS;
535: A->ops->setoption = MatSetOption_IS;
537: PetscMapSetBlockSize(&A->rmap,1);
538: PetscMapSetBlockSize(&A->cmap,1);
539: PetscMapSetUp(&A->rmap);
540: PetscMapSetUp(&A->cmap);
542: b->A = 0;
543: b->ctx = 0;
544: b->x = 0;
545: b->y = 0;
546: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);
547: PetscObjectChangeTypeName((PetscObject)A,MATIS);
549: return(0);
550: }