Actual source code: matis.c
petsc-3.3-p5 2012-12-01
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,bs;
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: MatGetBlockSize(A,&bs);
166: MatCreate(PETSC_COMM_SELF,&is->A);
167: MatSetSizes(is->A,n,n,n,n);
168: MatSetBlockSize(is->A,bs);
169: MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);
170: MatAppendOptionsPrefix(is->A,"is_");
171: MatSetFromOptions(is->A);
173: /* Create the local work vectors */
174: VecCreate(PETSC_COMM_SELF,&is->x);
175: VecSetBlockSize(is->x,bs);
176: VecSetSizes(is->x,n,n);
177: VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);
178: VecAppendOptionsPrefix(is->x,"is_");
179: VecSetFromOptions(is->x);
180: VecDuplicate(is->x,&is->y);
182: /* setup the global to local scatter */
183: ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);
184: ISLocalToGlobalMappingApplyIS(rmapping,to,&from);
185: MatGetVecs(A,&global,PETSC_NULL);
186: VecScatterCreate(global,from,is->x,to,&is->ctx);
187: VecDestroy(&global);
188: ISDestroy(&to);
189: ISDestroy(&from);
190: return(0);
191: }
193: #define ISG2LMapApply(mapping,n,in,out) 0;\
194: if (!(mapping)->globals) {\
195: PetscErrorCode _ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\
196: }\
197: {\
198: PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\
199: for (_i=0; _i<n; _i++) {\
200: if (in[_i] < 0) out[_i] = in[_i];\
201: else if (in[_i] < _start) out[_i] = -1;\
202: else if (in[_i] > _end) out[_i] = -1;\
203: else out[_i] = _globals[in[_i] - _start];\
204: }\
205: }
209: PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
210: {
211: Mat_IS *is = (Mat_IS*)mat->data;
212: PetscInt rows_l[2048],cols_l[2048];
216: #if defined(PETSC_USE_DEBUG)
217: if (m > 2048 || n > 2048) {
218: SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
219: }
220: #endif
221: ISG2LMapApply(is->mapping,m,rows,rows_l);
222: ISG2LMapApply(is->mapping,n,cols,cols_l);
223: MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);
224: return(0);
225: }
227: #undef ISG2LMapSetUp
228: #undef ISG2LMapApply
232: PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
233: {
235: Mat_IS *is = (Mat_IS*)A->data;
238: MatSetValues(is->A,m,rows,n,cols,values,addv);
239: return(0);
240: }
244: PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
245: {
246: Mat_IS *is = (Mat_IS*)A->data;
247: PetscInt n_l=0, *rows_l = PETSC_NULL;
251: if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
252: if (n) {
253: PetscMalloc(n*sizeof(PetscInt),&rows_l);
254: ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);
255: }
256: MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);
257: PetscFree(rows_l);
258: return(0);
259: }
263: PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
264: {
265: Mat_IS *is = (Mat_IS*)A->data;
267: PetscInt i;
268: PetscScalar *array;
271: if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
272: {
273: /*
274: Set up is->x as a "counting vector". This is in order to MatMult_IS
275: work properly in the interface nodes.
276: */
277: Vec counter;
278: PetscScalar one=1.0, zero=0.0;
279: MatGetVecs(A,&counter,PETSC_NULL);
280: VecSet(counter,zero);
281: VecSet(is->x,one);
282: VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);
283: VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);
284: VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);
285: VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);
286: VecDestroy(&counter);
287: }
288: if (!n) {
289: is->pure_neumann = PETSC_TRUE;
290: } else {
291: is->pure_neumann = PETSC_FALSE;
292: VecGetArray(is->x,&array);
293: MatZeroRows(is->A,n,rows,diag,0,0);
294: for (i=0; i<n; i++) {
295: MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);
296: }
297: MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);
298: MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);
299: VecRestoreArray(is->x,&array);
300: }
302: return(0);
303: }
307: PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
308: {
309: Mat_IS *is = (Mat_IS*)A->data;
313: MatAssemblyBegin(is->A,type);
314: return(0);
315: }
319: PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
320: {
321: Mat_IS *is = (Mat_IS*)A->data;
325: MatAssemblyEnd(is->A,type);
326: return(0);
327: }
329: EXTERN_C_BEGIN
332: PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
333: {
334: Mat_IS *is = (Mat_IS *)mat->data;
335:
337: *local = is->A;
338: return(0);
339: }
340: EXTERN_C_END
344: /*@
345: MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
347: Input Parameter:
348: . mat - the matrix
350: Output Parameter:
351: . local - the local matrix usually MATSEQAIJ
353: Level: advanced
355: Notes:
356: This can be called if you have precomputed the nonzero structure of the
357: matrix and want to provide it to the inner matrix object to improve the performance
358: of the MatSetValues() operation.
360: .seealso: MATIS
361: @*/
362: PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
363: {
369: PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));
370: return(0);
371: }
373: EXTERN_C_BEGIN
376: PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
377: {
378: Mat_IS *is = (Mat_IS *)mat->data;
379: PetscInt nrows,ncols,orows,ocols;
383: if(is->A) {
384: MatGetSize(is->A,&orows,&ocols);
385: MatGetSize(local,&nrows,&ncols);
386: if(orows != nrows || ocols != ncols ) {
387: SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols);
388: }
389: }
390: PetscObjectReference((PetscObject)local);
391: MatDestroy(&is->A);
392: is->A = local;
394: return(0);
395: }
396: EXTERN_C_END
400: /*@
401: MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
403: Input Parameter:
404: . mat - the matrix
405: . local - the local matrix usually MATSEQAIJ
407: Output Parameter:
409: Level: advanced
411: Notes:
412: This can be called if you have precomputed the local matrix and
413: want to provide it to the matrix object MATIS.
415: .seealso: MATIS
416: @*/
417: PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
418: {
423: PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));
424: return(0);
425: }
429: PetscErrorCode MatZeroEntries_IS(Mat A)
430: {
431: Mat_IS *a = (Mat_IS*)A->data;
435: MatZeroEntries(a->A);
436: return(0);
437: }
441: PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
442: {
443: Mat_IS *is = (Mat_IS*)A->data;
447: MatScale(is->A,a);
448: return(0);
449: }
453: PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
454: {
455: Mat_IS *is = (Mat_IS*)A->data;
459: /* get diagonal of the local matrix */
460: MatGetDiagonal(is->A,is->x);
462: /* scatter diagonal back into global vector */
463: VecSet(v,0);
464: VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);
465: VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);
466: return(0);
467: }
471: PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
472: {
473: Mat_IS *a = (Mat_IS*)A->data;
477: MatSetOption(a->A,op,flg);
478: return(0);
479: }
483: /*@
484: MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
485: process but not across processes.
487: Input Parameters:
488: + comm - MPI communicator that will share the matrix
489: . bs - local and global block size of the matrix
490: . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
491: - map - mapping that defines the global number for each local number
493: Output Parameter:
494: . A - the resulting matrix
496: Level: advanced
498: Notes: See MATIS for more details
499: m and n are NOT related to the size of the map, they are the size of the part of the vector owned
500: by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
501: plus the ghost points to global indices.
503: .seealso: MATIS, MatSetLocalToGlobalMapping()
504: @*/
505: PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
506: {
510: MatCreate(comm,A);
511: MatSetBlockSize(*A,bs);
512: MatSetSizes(*A,m,n,M,N);
513: MatSetType(*A,MATIS);
514: MatSetUp(*A);
515: MatSetLocalToGlobalMapping(*A,map,map);
516: return(0);
517: }
519: /*MC
520: MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
521: This stores the matrices in globally unassembled form. Each processor
522: assembles only its local Neumann problem and the parallel matrix vector
523: product is handled "implicitly".
525: Operations Provided:
526: + MatMult()
527: . MatMultAdd()
528: . MatMultTranspose()
529: . MatMultTransposeAdd()
530: . MatZeroEntries()
531: . MatSetOption()
532: . MatZeroRows()
533: . MatZeroRowsLocal()
534: . MatSetValues()
535: . MatSetValuesLocal()
536: . MatScale()
537: . MatGetDiagonal()
538: - MatSetLocalToGlobalMapping()
540: Options Database Keys:
541: . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
543: Notes: Options prefix for the inner matrix are given by -is_mat_xxx
544:
545: You must call MatSetLocalToGlobalMapping() before using this matrix type.
547: You can do matrix preallocation on the local matrix after you obtain it with
548: MatISGetLocalMat()
550: Level: advanced
552: .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
554: M*/
556: EXTERN_C_BEGIN
559: PetscErrorCode MatCreate_IS(Mat A)
560: {
562: Mat_IS *b;
565: PetscNewLog(A,Mat_IS,&b);
566: A->data = (void*)b;
567: PetscMemzero(A->ops,sizeof(struct _MatOps));
569: A->ops->mult = MatMult_IS;
570: A->ops->multadd = MatMultAdd_IS;
571: A->ops->multtranspose = MatMultTranspose_IS;
572: A->ops->multtransposeadd = MatMultTransposeAdd_IS;
573: A->ops->destroy = MatDestroy_IS;
574: A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
575: A->ops->setvalues = MatSetValues_IS;
576: A->ops->setvalueslocal = MatSetValuesLocal_IS;
577: A->ops->zerorows = MatZeroRows_IS;
578: A->ops->zerorowslocal = MatZeroRowsLocal_IS;
579: A->ops->assemblybegin = MatAssemblyBegin_IS;
580: A->ops->assemblyend = MatAssemblyEnd_IS;
581: A->ops->view = MatView_IS;
582: A->ops->zeroentries = MatZeroEntries_IS;
583: A->ops->scale = MatScale_IS;
584: A->ops->getdiagonal = MatGetDiagonal_IS;
585: A->ops->setoption = MatSetOption_IS;
587: PetscLayoutSetUp(A->rmap);
588: PetscLayoutSetUp(A->cmap);
590: b->A = 0;
591: b->ctx = 0;
592: b->x = 0;
593: b->y = 0;
594: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);
595: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);
596: PetscObjectChangeTypeName((PetscObject)A,MATIS);
598: return(0);
599: }
600: EXTERN_C_END