Actual source code: is.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/is.h

 17: int MatDestroy_IS(Mat A)
 18: {
 19:   int    ierr;
 20:   Mat_IS *b = (Mat_IS*)A->data;

 23:   if (b->A) {
 24:     MatDestroy(b->A);
 25:   }
 26:   if (b->ctx) {
 27:     VecScatterDestroy(b->ctx);
 28:   }
 29:   if (b->x) {
 30:     VecDestroy(b->x);
 31:   }
 32:   if (b->y) {
 33:     VecDestroy(b->y);
 34:   }
 35:   if (b->mapping) {
 36:     ISLocalToGlobalMappingDestroy(b->mapping);
 37:   }
 38:   PetscFree(b);
 39:   return(0);
 40: }

 42: int MatMult_IS(Mat A,Vec x,Vec y)
 43: {
 44:   int    ierr;
 45:   Mat_IS *is = (Mat_IS*)A->data;
 46:   PetscScalar zero = 0.0;

 49:   /*  scatter the global vector x into the local work vector */
 50:   VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
 51:   VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);

 53:   /* multiply the local matrix */
 54:   MatMult(is->A,is->x,is->y);

 56:   /* scatter product back into global memory */
 57:   VecSet(&zero,y);
 58:   VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);
 59:   VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);

 61:   return(0);
 62: }

 64: int MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
 65: {
 66:   int    ierr,n;
 67:   Mat_IS *is = (Mat_IS*)A->data;
 68:   IS     from,to;
 69:   Vec    global;

 72:   is->mapping = mapping;
 73:   PetscObjectReference((PetscObject)mapping);

 75:   /* Create the local matrix A */
 76:   ISLocalToGlobalMappingGetSize(mapping,&n);
 77:   MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);
 78:   PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");
 79:   MatSetFromOptions(is->A);

 81:   /* Create the local work vectors */
 82:   VecCreateSeq(PETSC_COMM_SELF,n,&is->x);
 83:   VecDuplicate(is->x,&is->y);

 85:   /* setup the global to local scatter */
 86:   ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);
 87:   ISLocalToGlobalMappingApplyIS(mapping,to,&from);
 88:   VecCreateMPI(A->comm,A->n,A->N,&global);
 89:   VecScatterCreate(global,from,is->x,to,&is->ctx);
 90:   VecDestroy(global);
 91:   ISDestroy(to);
 92:   ISDestroy(from);
 93:   return(0);
 94: }


 97: int MatSetValuesLocal_IS(Mat A,int m,int *rows,int n,int *cols,PetscScalar *values,InsertMode addv)
 98: {
 99:   int    ierr;
100:   Mat_IS *is = (Mat_IS*)A->data;

103:   MatSetValues(is->A,m,rows,n,cols,values,addv);
104:   return(0);
105: }

107: int MatZeroRowsLocal_IS(Mat A,IS isrows,PetscScalar *diag)
108: {
109:   Mat_IS *is = (Mat_IS*)A->data;
110:   int    ierr,i,n,*rows;
111:   PetscScalar *array;


115:   {
116:     /*
117:        Set up is->x as a "counting vector". This is in order to MatMult_IS
118:        work properly in the interface nodes.
119:     */
120:     Vec    counter;
121:     PetscScalar one=1.0, zero=0.0;
122:     VecCreateMPI(A->comm,A->n,A->N,&counter);
123:     VecSet(&zero,counter);
124:     VecSet(&one,is->x);
125:     VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);
126:     VecScatterEnd  (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);
127:     VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
128:     VecScatterEnd  (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);
129:     VecDestroy(counter);
130:   }
131:   ISGetLocalSize(isrows,&n);
132:   if (n == 0) {
133:     is->pure_neumann = PETSC_TRUE;
134:   } else {
135:     is->pure_neumann = PETSC_FALSE;
136:     ISGetIndices(isrows,&rows);
137:     VecGetArray(is->x,&array);
138:     MatZeroRows(is->A,isrows,diag);
139:     for (i=0; i<n; i++) {
140:       MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);
141:     }
142:     MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);
143:     MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);
144:     VecRestoreArray(is->x,&array);
145:     ISRestoreIndices(isrows,&rows);
146:   }

148:   return(0);
149: }

151: int MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
152: {
153:   Mat_IS *is = (Mat_IS*)A->data;
154:   int    ierr;
156:   MatAssemblyBegin(is->A,type);
157:   return(0);
158: }

160: int MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
161: {
162:   Mat_IS *is = (Mat_IS*)A->data;
163:   int    ierr;
165:   MatAssemblyEnd(is->A,type);
166:   return(0);
167: }

169: EXTERN_C_BEGIN
170: int MatCreate_IS(Mat A)
171: {
172:   int    ierr;
173:   Mat_IS *b;

176:   ierr                = PetscNew(Mat_IS,&b);
177:   A->data             = (void*)b;
178:   PetscMemzero(b,sizeof(Mat_IS));
179:   PetscMemzero(A->ops,sizeof(struct _MatOps));
180:   A->factor           = 0;
181:   A->mapping          = 0;

183:   A->ops->mult                    = MatMult_IS;
184:   A->ops->destroy                 = MatDestroy_IS;
185:   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
186:   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
187:   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
188:   A->ops->assemblybegin           = MatAssemblyBegin_IS;
189:   A->ops->assemblyend             = MatAssemblyEnd_IS;

191:   PetscSplitOwnership(A->comm,&A->m,&A->M);
192:   PetscSplitOwnership(A->comm,&A->n,&A->N);
193:   MPI_Scan(&A->m,&b->rend,1,MPI_INT,MPI_SUM,A->comm);
194:   b->rstart = b->rend - A->m;

196:   b->A          = 0;
197:   b->ctx        = 0;
198:   b->x          = 0;
199:   b->y          = 0;

201:   return(0);
202: }
203: EXTERN_C_END