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