1 /* 2 Creates a matrix class for using the Neumann-Neumann type preconditioners. 3 This stores the matrices in globally unassembled form. Each processor 4 assembles only its local Neumann problem and the parallel matrix vector 5 product is handled "implicitly". 6 7 We provide: 8 MatMult() 9 10 Currently this allows for only one subdomain per processor. 11 12 */ 13 14 #include "src/mat/impls/is/matis.h" /*I "petscmat.h" I*/ 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "MatDestroy_IS" 18 int MatDestroy_IS(Mat A) 19 { 20 int ierr; 21 Mat_IS *b = (Mat_IS*)A->data; 22 23 PetscFunctionBegin; 24 if (b->A) { 25 ierr = MatDestroy(b->A);CHKERRQ(ierr); 26 } 27 if (b->ctx) { 28 ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 29 } 30 if (b->x) { 31 ierr = VecDestroy(b->x);CHKERRQ(ierr); 32 } 33 if (b->y) { 34 ierr = VecDestroy(b->y);CHKERRQ(ierr); 35 } 36 if (b->mapping) { 37 ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr); 38 } 39 ierr = PetscFree(b);CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "MatMult_IS" 45 int MatMult_IS(Mat A,Vec x,Vec y) 46 { 47 int ierr; 48 Mat_IS *is = (Mat_IS*)A->data; 49 PetscScalar zero = 0.0; 50 51 PetscFunctionBegin; 52 /* scatter the global vector x into the local work vector */ 53 ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 54 ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 55 56 /* multiply the local matrix */ 57 ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 58 59 /* scatter product back into global memory */ 60 ierr = VecSet(&zero,y);CHKERRQ(ierr); 61 ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 62 ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 63 64 PetscFunctionReturn(0); 65 } 66 67 #undef __FUNCT__ 68 #define __FUNCT__ "MatView_IS" 69 int MatView_IS(Mat A,PetscViewer viewer) 70 { 71 Mat_IS *a = (Mat_IS*)A->data; 72 int ierr; 73 PetscViewer sviewer; 74 75 PetscFunctionBegin; 76 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 77 ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 78 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 84 int MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping) 85 { 86 int ierr,n; 87 Mat_IS *is = (Mat_IS*)A->data; 88 IS from,to; 89 Vec global; 90 91 PetscFunctionBegin; 92 is->mapping = mapping; 93 ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); 94 95 /* Create the local matrix A */ 96 ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 97 ierr = MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);CHKERRQ(ierr); 98 ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr); 99 ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 100 101 /* Create the local work vectors */ 102 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr); 103 ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 104 105 /* setup the global to local scatter */ 106 ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 107 ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr); 108 ierr = VecCreateMPI(A->comm,A->n,A->N,&global);CHKERRQ(ierr); 109 ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 110 ierr = VecDestroy(global);CHKERRQ(ierr); 111 ierr = ISDestroy(to);CHKERRQ(ierr); 112 ierr = ISDestroy(from);CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 117 #undef __FUNCT__ 118 #define __FUNCT__ "MatSetValuesLocal_IS" 119 int MatSetValuesLocal_IS(Mat A,int m,const int *rows, int n,const int *cols,const PetscScalar *values,InsertMode addv) 120 { 121 int ierr; 122 Mat_IS *is = (Mat_IS*)A->data; 123 124 PetscFunctionBegin; 125 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "MatZeroRowsLocal_IS" 131 int MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag) 132 { 133 Mat_IS *is = (Mat_IS*)A->data; 134 int ierr,i,n,*rows; 135 PetscScalar *array; 136 137 PetscFunctionBegin; 138 139 { 140 /* 141 Set up is->x as a "counting vector". This is in order to MatMult_IS 142 work properly in the interface nodes. 143 */ 144 Vec counter; 145 PetscScalar one=1.0, zero=0.0; 146 ierr = VecCreateMPI(A->comm,A->n,A->N,&counter);CHKERRQ(ierr); 147 ierr = VecSet(&zero,counter);CHKERRQ(ierr); 148 ierr = VecSet(&one,is->x);CHKERRQ(ierr); 149 ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 150 ierr = VecScatterEnd (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 151 ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 152 ierr = VecScatterEnd (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 153 ierr = VecDestroy(counter);CHKERRQ(ierr); 154 } 155 ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr); 156 if (!n) { 157 is->pure_neumann = PETSC_TRUE; 158 } else { 159 is->pure_neumann = PETSC_FALSE; 160 ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 161 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 162 ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr); 163 for (i=0; i<n; i++) { 164 ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 165 } 166 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 167 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 169 ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 170 } 171 172 PetscFunctionReturn(0); 173 } 174 175 #undef __FUNCT__ 176 #define __FUNCT__ "MatAssemblyBegin_IS" 177 int MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 178 { 179 Mat_IS *is = (Mat_IS*)A->data; 180 int ierr; 181 182 PetscFunctionBegin; 183 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 #undef __FUNCT__ 188 #define __FUNCT__ "MatAssemblyEnd_IS" 189 int MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 190 { 191 Mat_IS *is = (Mat_IS*)A->data; 192 int ierr; 193 194 PetscFunctionBegin; 195 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 196 PetscFunctionReturn(0); 197 } 198 199 EXTERN_C_BEGIN 200 #undef __FUNCT__ 201 #define __FUNCT__ "MatISGetLocalMat_IS" 202 int MatISGetLocalMat_IS(Mat mat,Mat *local) 203 { 204 Mat_IS *is = (Mat_IS *)mat->data; 205 206 PetscFunctionBegin; 207 *local = is->A; 208 PetscFunctionReturn(0); 209 } 210 EXTERN_C_END 211 212 #undef __FUNCT__ 213 #define __FUNCT__ "MatISGetLocalMat" 214 /*@ 215 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 216 217 Input Parameter: 218 . mat - the matrix 219 220 Output Parameter: 221 . local - the local matrix usually MATSEQAIJ 222 223 Level: advanced 224 225 Notes: 226 This can be called if you have precomputed the nonzero structure of the 227 matrix and want to provide it to the inner matrix object to improve the performance 228 of the MatSetValues() operation. 229 230 .seealso: MATIS 231 @*/ 232 int MatISGetLocalMat(Mat mat,Mat *local) 233 { 234 int ierr,(*f)(Mat,Mat *); 235 236 PetscFunctionBegin; 237 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 238 PetscValidPointer(local,2); 239 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr); 240 if (f) { 241 ierr = (*f)(mat,local);CHKERRQ(ierr); 242 } else { 243 local = 0; 244 } 245 PetscFunctionReturn(0); 246 } 247 248 /*MC 249 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 250 This stores the matrices in globally unassembled form. Each processor 251 assembles only its local Neumann problem and the parallel matrix vector 252 product is handled "implicitly". 253 254 Operations Provided: 255 . MatMult 256 257 Options Database Keys: 258 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 259 260 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 261 262 You must call MatSetLocalToGlobalMapping() before using this matrix type. 263 264 You can do matrix preallocation on the local matrix after you obtain it with 265 MatISGetLocalMat() 266 267 Level: advanced 268 269 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 270 271 M*/ 272 273 EXTERN_C_BEGIN 274 #undef __FUNCT__ 275 #define __FUNCT__ "MatCreate_IS" 276 int MatCreate_IS(Mat A) 277 { 278 int ierr; 279 Mat_IS *b; 280 281 PetscFunctionBegin; 282 ierr = PetscNew(Mat_IS,&b);CHKERRQ(ierr); 283 A->data = (void*)b; 284 ierr = PetscMemzero(b,sizeof(Mat_IS));CHKERRQ(ierr); 285 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 286 A->factor = 0; 287 A->mapping = 0; 288 289 A->ops->mult = MatMult_IS; 290 A->ops->destroy = MatDestroy_IS; 291 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 292 A->ops->setvalueslocal = MatSetValuesLocal_IS; 293 A->ops->zerorowslocal = MatZeroRowsLocal_IS; 294 A->ops->assemblybegin = MatAssemblyBegin_IS; 295 A->ops->assemblyend = MatAssemblyEnd_IS; 296 A->ops->view = MatView_IS; 297 298 ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr); 299 ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr); 300 ierr = MPI_Scan(&A->m,&b->rend,1,MPI_INT,MPI_SUM,A->comm);CHKERRQ(ierr); 301 b->rstart = b->rend - A->m; 302 303 b->A = 0; 304 b->ctx = 0; 305 b->x = 0; 306 b->y = 0; 307 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 308 309 PetscFunctionReturn(0); 310 } 311 EXTERN_C_END 312 313 314 315 316 317 318 319 320 321 322 323 324 325