Actual source code: matis.c
  1: /*$Id: is.c,v 1.15 2001/08/06 21:15:46 bsmith Exp $*/
  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
 19: int MatDestroy_IS(Mat A)
 20: {
 21:   int    ierr;
 22:   Mat_IS *b = (Mat_IS*)A->data;
 25:   if (b->A) {
 26:     MatDestroy(b->A);
 27:   }
 28:   if (b->ctx) {
 29:     VecScatterDestroy(b->ctx);
 30:   }
 31:   if (b->x) {
 32:     VecDestroy(b->x);
 33:   }
 34:   if (b->y) {
 35:     VecDestroy(b->y);
 36:   }
 37:   if (b->mapping) {
 38:     ISLocalToGlobalMappingDestroy(b->mapping);
 39:   }
 40:   PetscFree(b);
 41:   return(0);
 42: }
 46: int MatMult_IS(Mat A,Vec x,Vec y)
 47: {
 48:   int         ierr;
 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: int MatView_IS(Mat A,PetscViewer viewer)
 71: {
 72:   Mat_IS      *a = (Mat_IS*)A->data;
 73:   int         ierr;
 74:   PetscViewer sviewer;
 77:   PetscViewerGetSingleton(viewer,&sviewer);
 78:   MatView(a->A,sviewer);
 79:   PetscViewerRestoreSingleton(viewer,&sviewer);
 80:   return(0);
 81: }
 85: int MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
 86: {
 87:   int    ierr,n;
 88:   Mat_IS *is = (Mat_IS*)A->data;
 89:   IS     from,to;
 90:   Vec    global;
 93:   is->mapping = mapping;
 94:   PetscObjectReference((PetscObject)mapping);
 96:   /* Create the local matrix A */
 97:   ISLocalToGlobalMappingGetSize(mapping,&n);
 98:   MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);
 99:   PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");
100:   MatSetFromOptions(is->A);
102:   /* Create the local work vectors */
103:   VecCreateSeq(PETSC_COMM_SELF,n,&is->x);
104:   VecDuplicate(is->x,&is->y);
106:   /* setup the global to local scatter */
107:   ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);
108:   ISLocalToGlobalMappingApplyIS(mapping,to,&from);
109:   VecCreateMPI(A->comm,A->n,A->N,&global);
110:   VecScatterCreate(global,from,is->x,to,&is->ctx);
111:   VecDestroy(global);
112:   ISDestroy(to);
113:   ISDestroy(from);
114:   return(0);
115: }
120: int MatSetValuesLocal_IS(Mat A,int m,const int *rows, int n,const int *cols,const PetscScalar *values,InsertMode addv)
121: {
122:   int    ierr;
123:   Mat_IS *is = (Mat_IS*)A->data;
126:   MatSetValues(is->A,m,rows,n,cols,values,addv);
127:   return(0);
128: }
132: int MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag)
133: {
134:   Mat_IS      *is = (Mat_IS*)A->data;
135:   int         ierr,i,n,*rows;
136:   PetscScalar *array;
140:   {
141:     /*
142:        Set up is->x as a "counting vector". This is in order to MatMult_IS
143:        work properly in the interface nodes.
144:     */
145:     Vec         counter;
146:     PetscScalar one=1.0, zero=0.0;
147:     VecCreateMPI(A->comm,A->n,A->N,&counter);
148:     VecSet(&zero,counter);
149:     VecSet(&one,is->x);
150:     VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);
151:     VecScatterEnd  (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);
152:     VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
153:     VecScatterEnd  (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
154:     VecDestroy(counter);
155:   }
156:   ISGetLocalSize(isrows,&n);
157:   if (n == 0) {
158:     is->pure_neumann = PETSC_TRUE;
159:   } else {
160:     is->pure_neumann = PETSC_FALSE;
161:     ISGetIndices(isrows,&rows);
162:     VecGetArray(is->x,&array);
163:     MatZeroRows(is->A,isrows,diag);
164:     for (i=0; i<n; i++) {
165:       MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);
166:     }
167:     MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);
168:     MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);
169:     VecRestoreArray(is->x,&array);
170:     ISRestoreIndices(isrows,&rows);
171:   }
173:   return(0);
174: }
178: int MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
179: {
180:   Mat_IS *is = (Mat_IS*)A->data;
181:   int    ierr;
184:   MatAssemblyBegin(is->A,type);
185:   return(0);
186: }
190: int MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
191: {
192:   Mat_IS *is = (Mat_IS*)A->data;
193:   int    ierr;
196:   MatAssemblyEnd(is->A,type);
197:   return(0);
198: }
200: EXTERN_C_BEGIN
203: int MatISGetLocalMat_IS(Mat mat,Mat *local)
204: {
205:   Mat_IS *is = (Mat_IS *)mat->data;
206: 
208:   *local = is->A;
209:   return(0);
210: }
211: EXTERN_C_END
215: /*@
216:     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
218:   Input Parameter:
219: .  mat - the matrix
221:   Output Parameter:
222: .  local - the local matrix usually MATSEQAIJ
224:   Level: advanced
226:   Notes:
227:     This can be called if you have precomputed the nonzero structure of the 
228:   matrix and want to provide it to the inner matrix object to improve the performance
229:   of the MatSetValues() operation.
231: .seealso: MATIS
232: @*/
233: int MatISGetLocalMat(Mat mat,Mat *local)
234: {
235:   int ierr,(*f)(Mat,Mat *);
240:   PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);
241:   if (f) {
242:     (*f)(mat,local);
243:   } else {
244:     local = 0;
245:   }
246:   return(0);
247: }
249: /*MC
250:    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
251:    This stores the matrices in globally unassembled form. Each processor 
252:    assembles only its local Neumann problem and the parallel matrix vector 
253:    product is handled "implicitly".
255:    Operations Provided:
256: .  MatMult
258:    Options Database Keys:
259: . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
261:    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
262:     
263:           You must call MatSetLocalToGlobalMapping() before using this matrix type.
265:           You can do matrix preallocation on the local matrix after you obtain it with 
266:           MatISGetLocalMat()
268:   Level: advanced
270: .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
272: M*/
274: EXTERN_C_BEGIN
277: int MatCreate_IS(Mat A)
278: {
279:   int    ierr;
280:   Mat_IS *b;
283:   PetscNew(Mat_IS,&b);
284:   A->data             = (void*)b;
285:   PetscMemzero(b,sizeof(Mat_IS));
286:   PetscMemzero(A->ops,sizeof(struct _MatOps));
287:   A->factor           = 0;
288:   A->mapping          = 0;
290:   A->ops->mult                    = MatMult_IS;
291:   A->ops->destroy                 = MatDestroy_IS;
292:   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
293:   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
294:   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
295:   A->ops->assemblybegin           = MatAssemblyBegin_IS;
296:   A->ops->assemblyend             = MatAssemblyEnd_IS;
297:   A->ops->view                    = MatView_IS;
299:   PetscSplitOwnership(A->comm,&A->m,&A->M);
300:   PetscSplitOwnership(A->comm,&A->n,&A->N);
301:   MPI_Scan(&A->m,&b->rend,1,MPI_INT,MPI_SUM,A->comm);
302:   b->rstart = b->rend - A->m;
304:   b->A          = 0;
305:   b->ctx        = 0;
306:   b->x          = 0;
307:   b->y          = 0;
308:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);
310:   return(0);
311: }
312: EXTERN_C_END