Actual source code: matis.c
2: /*
3: Creates a matrix class for using the Neumann-Neumann type preconditioners.
4: This stores the matrices in globally unassembled form. Each processor
5: assembles only its local Neumann problem and the parallel matrix vector
6: product is handled "implicitly".
8: We provide:
9: MatMult()
11: Currently this allows for only one subdomain per processor.
13: */
15: #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/
19: PetscErrorCode MatDestroy_IS(Mat A)
20: {
22: Mat_IS *b = (Mat_IS*)A->data;
25: MatDestroy(&b->A);
26: VecScatterDestroy(&b->ctx);
27: VecDestroy(&b->x);
28: VecDestroy(&b->y);
29: ISLocalToGlobalMappingDestroy(&b->mapping);
30: PetscFree(A->data);
31: PetscObjectChangeTypeName((PetscObject)A,0);
32: PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);
33: return(0);
34: }
38: PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
39: {
41: Mat_IS *is = (Mat_IS*)A->data;
42: PetscScalar zero = 0.0;
45: /* scatter the global vector x into the local work vector */
46: VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
47: VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
49: /* multiply the local matrix */
50: MatMult(is->A,is->x,is->y);
52: /* scatter product back into global memory */
53: VecSet(y,zero);
54: VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
55: VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
57: return(0);
58: }
62: static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
63: {
65: Mat_IS *is = (Mat_IS*)A->data;
68: /* scatter the global vector v1 into the local work vector */
69: VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
70: VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
71: VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
72: VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
74: /* multiply the local matrix */
75: MatMultAdd(is->A,is->x,is->y,is->y);
77: /* scatter result back into global vector */
78: VecSet(v3,0);
79: VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
80: VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
81: return(0);
82: }
86: PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
87: {
88: Mat_IS *is = (Mat_IS*)A->data;
92: /* scatter the global vector x into the local work vector */
93: VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
94: VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);
96: /* multiply the local matrix */
97: MatMultTranspose(is->A,is->x,is->y);
99: /* scatter product back into global vector */
100: VecSet(y,0);
101: VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
102: VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);
103: return(0);
104: }
108: PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
109: {
110: Mat_IS *is = (Mat_IS*)A->data;
114: /* scatter the global vectors v1 and v2 into the local work vectors */
115: VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
116: VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);
117: VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
118: VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);
120: /* multiply the local matrix */
121: MatMultTransposeAdd(is->A,is->x,is->y,is->y);
123: /* scatter result back into global vector */
124: VecSet(v3,0);
125: VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
126: VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);
127: return(0);
128: }
132: PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
133: {
134: Mat_IS *a = (Mat_IS*)A->data;
136: PetscViewer sviewer;
139: PetscViewerGetSingleton(viewer,&sviewer);
140: MatView(a->A,sviewer);
141: PetscViewerRestoreSingleton(viewer,&sviewer);
142: return(0);
143: }
147: PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
148: {
150: PetscInt n;
151: Mat_IS *is = (Mat_IS*)A->data;
152: IS from,to;
153: Vec global;
156: if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
158: if (rmapping != cmapping) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
159: PetscObjectReference((PetscObject)rmapping);
160: ISLocalToGlobalMappingDestroy(&is->mapping);
161: is->mapping = rmapping;
163: /* Create the local matrix A */
164: ISLocalToGlobalMappingGetSize(rmapping,&n);
165: MatCreate(PETSC_COMM_SELF,&is->A);
166: MatSetSizes(is->A,n,n,n,n);
167: MatSetOptionsPrefix(is->A,"is");
168: MatSetFromOptions(is->A);
170: /* Create the local work vectors */
171: VecCreateSeq(PETSC_COMM_SELF,n,&is->x);
172: VecDuplicate(is->x,&is->y);
174: /* setup the global to local scatter */
175: ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);
176: ISLocalToGlobalMappingApplyIS(rmapping,to,&from);
177: VecCreateMPIWithArray(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,PETSC_NULL,&global);
178: VecScatterCreate(global,from,is->x,to,&is->ctx);
179: VecDestroy(&global);
180: ISDestroy(&to);
181: ISDestroy(&from);
182: return(0);
183: }
185: #define ISG2LMapApply(mapping,n,in,out) 0;\
186: if (!(mapping)->globals) {\
187: PetscErrorCode _ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\
188: }\
189: {\
190: PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\
191: for (_i=0; _i<n; _i++) {\
192: if (in[_i] < 0) out[_i] = in[_i];\
193: else if (in[_i] < _start) out[_i] = -1;\
194: else if (in[_i] > _end) out[_i] = -1;\
195: else out[_i] = _globals[in[_i] - _start];\
196: }\
197: }
201: PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
202: {
203: Mat_IS *is = (Mat_IS*)mat->data;
204: PetscInt rows_l[2048],cols_l[2048];
208: #if defined(PETSC_USE_DEBUG)
209: if (m > 2048 || n > 2048) {
210: SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
211: }
212: #endif
213: ISG2LMapApply(is->mapping,m,rows,rows_l);
214: ISG2LMapApply(is->mapping,n,cols,cols_l);
215: MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);
216: return(0);
217: }
219: #undef ISG2LMapSetUp
220: #undef ISG2LMapApply
224: PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
225: {
227: Mat_IS *is = (Mat_IS*)A->data;
230: MatSetValues(is->A,m,rows,n,cols,values,addv);
231: return(0);
232: }
236: PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
237: {
238: Mat_IS *is = (Mat_IS*)A->data;
239: PetscInt n_l=0, *rows_l = PETSC_NULL;
243: if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
244: if (n) {
245: PetscMalloc(n*sizeof(PetscInt),&rows_l);
246: ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);
247: }
248: MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);
249: PetscFree(rows_l);
250: return(0);
251: }
255: PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
256: {
257: Mat_IS *is = (Mat_IS*)A->data;
259: PetscInt i;
260: PetscScalar *array;
263: if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
264: {
265: /*
266: Set up is->x as a "counting vector". This is in order to MatMult_IS
267: work properly in the interface nodes.
268: */
269: Vec counter;
270: PetscScalar one=1.0, zero=0.0;
271: VecCreateMPI(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,&counter);
272: VecSet(counter,zero);
273: VecSet(is->x,one);
274: VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);
275: VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);
276: VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);
277: VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);
278: VecDestroy(&counter);
279: }
280: if (!n) {
281: is->pure_neumann = PETSC_TRUE;
282: } else {
283: is->pure_neumann = PETSC_FALSE;
284: VecGetArray(is->x,&array);
285: MatZeroRows(is->A,n,rows,diag,0,0);
286: for (i=0; i<n; i++) {
287: MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);
288: }
289: MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);
290: MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);
291: VecRestoreArray(is->x,&array);
292: }
294: return(0);
295: }
299: PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
300: {
301: Mat_IS *is = (Mat_IS*)A->data;
305: MatAssemblyBegin(is->A,type);
306: return(0);
307: }
311: PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
312: {
313: Mat_IS *is = (Mat_IS*)A->data;
317: MatAssemblyEnd(is->A,type);
318: return(0);
319: }
324: PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
325: {
326: Mat_IS *is = (Mat_IS *)mat->data;
327:
329: *local = is->A;
330: return(0);
331: }
336: /*@
337: MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
339: Input Parameter:
340: . mat - the matrix
342: Output Parameter:
343: . local - the local matrix usually MATSEQAIJ
345: Level: advanced
347: Notes:
348: This can be called if you have precomputed the nonzero structure of the
349: matrix and want to provide it to the inner matrix object to improve the performance
350: of the MatSetValues() operation.
352: .seealso: MATIS
353: @*/
354: PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
355: {
361: PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));
362: return(0);
363: }
367: PetscErrorCode MatZeroEntries_IS(Mat A)
368: {
369: Mat_IS *a = (Mat_IS*)A->data;
373: MatZeroEntries(a->A);
374: return(0);
375: }
379: PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
380: {
381: Mat_IS *is = (Mat_IS*)A->data;
385: MatScale(is->A,a);
386: return(0);
387: }
391: PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
392: {
393: Mat_IS *is = (Mat_IS*)A->data;
397: /* get diagonal of the local matrix */
398: MatGetDiagonal(is->A,is->x);
400: /* scatter diagonal back into global vector */
401: VecSet(v,0);
402: VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);
403: VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);
404: return(0);
405: }
409: PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
410: {
411: Mat_IS *a = (Mat_IS*)A->data;
415: MatSetOption(a->A,op,flg);
416: return(0);
417: }
421: /*@
422: MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
423: process but not across processes.
425: Input Parameters:
426: + comm - MPI communicator that will share the matrix
427: . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
428: - map - mapping that defines the global number for each local number
430: Output Parameter:
431: . A - the resulting matrix
433: Level: advanced
435: Notes: See MATIS for more details
436: m and n are NOT related to the size of the map, they are the size of the part of the vector owned
437: by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
438: plus the ghost points to global indices.
440: .seealso: MATIS, MatSetLocalToGlobalMapping()
441: @*/
442: PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
443: {
447: MatCreate(comm,A);
448: MatSetSizes(*A,m,n,M,N);
449: MatSetType(*A,MATIS);
450: MatSetLocalToGlobalMapping(*A,map,map);
451: return(0);
452: }
454: /*MC
455: MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
456: This stores the matrices in globally unassembled form. Each processor
457: assembles only its local Neumann problem and the parallel matrix vector
458: product is handled "implicitly".
460: Operations Provided:
461: + MatMult()
462: . MatMultAdd()
463: . MatMultTranspose()
464: . MatMultTransposeAdd()
465: . MatZeroEntries()
466: . MatSetOption()
467: . MatZeroRows()
468: . MatZeroRowsLocal()
469: . MatSetValues()
470: . MatSetValuesLocal()
471: . MatScale()
472: . MatGetDiagonal()
473: - MatSetLocalToGlobalMapping()
475: Options Database Keys:
476: . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
478: Notes: Options prefix for the inner matrix are given by -is_mat_xxx
479:
480: You must call MatSetLocalToGlobalMapping() before using this matrix type.
482: You can do matrix preallocation on the local matrix after you obtain it with
483: MatISGetLocalMat()
485: Level: advanced
487: .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
489: M*/
494: PetscErrorCode MatCreate_IS(Mat A)
495: {
497: Mat_IS *b;
500: PetscNewLog(A,Mat_IS,&b);
501: A->data = (void*)b;
502: PetscMemzero(A->ops,sizeof(struct _MatOps));
504: A->ops->mult = MatMult_IS;
505: A->ops->multadd = MatMultAdd_IS;
506: A->ops->multtranspose = MatMultTranspose_IS;
507: A->ops->multtransposeadd = MatMultTransposeAdd_IS;
508: A->ops->destroy = MatDestroy_IS;
509: A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
510: A->ops->setvalues = MatSetValues_IS;
511: A->ops->setvalueslocal = MatSetValuesLocal_IS;
512: A->ops->zerorows = MatZeroRows_IS;
513: A->ops->zerorowslocal = MatZeroRowsLocal_IS;
514: A->ops->assemblybegin = MatAssemblyBegin_IS;
515: A->ops->assemblyend = MatAssemblyEnd_IS;
516: A->ops->view = MatView_IS;
517: A->ops->zeroentries = MatZeroEntries_IS;
518: A->ops->scale = MatScale_IS;
519: A->ops->getdiagonal = MatGetDiagonal_IS;
520: A->ops->setoption = MatSetOption_IS;
522: PetscLayoutSetBlockSize(A->rmap,1);
523: PetscLayoutSetBlockSize(A->cmap,1);
524: PetscLayoutSetUp(A->rmap);
525: PetscLayoutSetUp(A->cmap);
527: b->A = 0;
528: b->ctx = 0;
529: b->x = 0;
530: b->y = 0;
531: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);
532: PetscObjectChangeTypeName((PetscObject)A,MATIS);
534: return(0);
535: }