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". 7 8 We provide: 9 MatMult() 10 11 Currently this allows for only one subdomain per processor. 12 13 */ 14 15 #include "src/mat/impls/is/matis.h" /*I "petscmat.h" I*/ 16 17 #undef __FUNCT__ 18 #define __FUNCT__ "MatDestroy_IS" 19 int MatDestroy_IS(Mat A) 20 { 21 int ierr; 22 Mat_IS *b = (Mat_IS*)A->data; 23 24 PetscFunctionBegin; 25 if (b->A) { 26 ierr = MatDestroy(b->A);CHKERRQ(ierr); 27 } 28 if (b->ctx) { 29 ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 30 } 31 if (b->x) { 32 ierr = VecDestroy(b->x);CHKERRQ(ierr); 33 } 34 if (b->y) { 35 ierr = VecDestroy(b->y);CHKERRQ(ierr); 36 } 37 if (b->mapping) { 38 ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr); 39 } 40 ierr = PetscFree(b);CHKERRQ(ierr); 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "MatMult_IS" 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; 51 52 PetscFunctionBegin; 53 /* scatter the global vector x into the local work vector */ 54 ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 55 ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 56 57 /* multiply the local matrix */ 58 ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 59 60 /* scatter product back into global memory */ 61 ierr = VecSet(&zero,y);CHKERRQ(ierr); 62 ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 63 ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 64 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "MatView_IS" 70 int MatView_IS(Mat A,PetscViewer viewer) 71 { 72 Mat_IS *a = (Mat_IS*)A->data; 73 int ierr; 74 PetscViewer sviewer; 75 76 PetscFunctionBegin; 77 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 78 ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 79 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 80 PetscFunctionReturn(0); 81 } 82 83 #undef __FUNCT__ 84 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 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; 91 92 PetscFunctionBegin; 93 is->mapping = mapping; 94 ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); 95 96 /* Create the local matrix A */ 97 ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 98 ierr = MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);CHKERRQ(ierr); 99 ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr); 100 ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 101 102 /* Create the local work vectors */ 103 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr); 104 ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 105 106 /* setup the global to local scatter */ 107 ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 108 ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr); 109 ierr = VecCreateMPI(A->comm,A->n,A->N,&global);CHKERRQ(ierr); 110 ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 111 ierr = VecDestroy(global);CHKERRQ(ierr); 112 ierr = ISDestroy(to);CHKERRQ(ierr); 113 ierr = ISDestroy(from);CHKERRQ(ierr); 114 PetscFunctionReturn(0); 115 } 116 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "MatSetValuesLocal_IS" 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; 124 125 PetscFunctionBegin; 126 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNCT__ 131 #define __FUNCT__ "MatZeroRowsLocal_IS" 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; 137 138 PetscFunctionBegin; 139 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 ierr = VecCreateMPI(A->comm,A->n,A->N,&counter);CHKERRQ(ierr); 148 ierr = VecSet(&zero,counter);CHKERRQ(ierr); 149 ierr = VecSet(&one,is->x);CHKERRQ(ierr); 150 ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 151 ierr = VecScatterEnd (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 152 ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 153 ierr = VecScatterEnd (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 154 ierr = VecDestroy(counter);CHKERRQ(ierr); 155 } 156 ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr); 157 if (n == 0) { 158 is->pure_neumann = PETSC_TRUE; 159 } else { 160 is->pure_neumann = PETSC_FALSE; 161 ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 162 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 163 ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr); 164 for (i=0; i<n; i++) { 165 ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 166 } 167 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 170 ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 171 } 172 173 PetscFunctionReturn(0); 174 } 175 176 #undef __FUNCT__ 177 #define __FUNCT__ "MatAssemblyBegin_IS" 178 int MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 179 { 180 Mat_IS *is = (Mat_IS*)A->data; 181 int ierr; 182 183 PetscFunctionBegin; 184 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 185 PetscFunctionReturn(0); 186 } 187 188 #undef __FUNCT__ 189 #define __FUNCT__ "MatAssemblyEnd_IS" 190 int MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 191 { 192 Mat_IS *is = (Mat_IS*)A->data; 193 int ierr; 194 195 PetscFunctionBegin; 196 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 197 PetscFunctionReturn(0); 198 } 199 200 EXTERN_C_BEGIN 201 #undef __FUNCT__ 202 #define __FUNCT__ "MatISGetLocalMat_IS" 203 int MatISGetLocalMat_IS(Mat mat,Mat *local) 204 { 205 Mat_IS *is = (Mat_IS *)mat->data; 206 207 PetscFunctionBegin; 208 *local = is->A; 209 PetscFunctionReturn(0); 210 } 211 EXTERN_C_END 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "MatISGetLocalMat" 215 /*@ 216 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 217 218 Input Parameter: 219 . mat - the matrix 220 221 Output Parameter: 222 . local - the local matrix usually MATSEQAIJ 223 224 Level: advanced 225 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. 230 231 .seealso: MATIS 232 @*/ 233 int MatISGetLocalMat(Mat mat,Mat *local) 234 { 235 int ierr,(*f)(Mat,Mat *); 236 237 PetscFunctionBegin; 238 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 239 PetscValidPointer(local,2); 240 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr); 241 if (f) { 242 ierr = (*f)(mat,local);CHKERRQ(ierr); 243 } else { 244 local = 0; 245 } 246 PetscFunctionReturn(0); 247 } 248 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". 254 255 Operations Provided: 256 . MatMult 257 258 Options Database Keys: 259 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 260 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. 264 265 You can do matrix preallocation on the local matrix after you obtain it with 266 MatISGetLocalMat() 267 268 Level: advanced 269 270 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 271 272 M*/ 273 274 EXTERN_C_BEGIN 275 #undef __FUNCT__ 276 #define __FUNCT__ "MatCreate_IS" 277 int MatCreate_IS(Mat A) 278 { 279 int ierr; 280 Mat_IS *b; 281 282 PetscFunctionBegin; 283 ierr = PetscNew(Mat_IS,&b);CHKERRQ(ierr); 284 A->data = (void*)b; 285 ierr = PetscMemzero(b,sizeof(Mat_IS));CHKERRQ(ierr); 286 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 287 A->factor = 0; 288 A->mapping = 0; 289 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; 298 299 ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr); 300 ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr); 301 ierr = MPI_Scan(&A->m,&b->rend,1,MPI_INT,MPI_SUM,A->comm);CHKERRQ(ierr); 302 b->rstart = b->rend - A->m; 303 304 b->A = 0; 305 b->ctx = 0; 306 b->x = 0; 307 b->y = 0; 308 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 309 310 PetscFunctionReturn(0); 311 } 312 EXTERN_C_END 313 314 315 316 317 318 319 320 321 322 323 324 325 326