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,IS isrows,const PetscScalar *diag) 137 { 138 Mat_IS *is = (Mat_IS*)A->data; 139 PetscErrorCode ierr; 140 PetscInt i,n,*rows; 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 ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr); 161 if (!n) { 162 is->pure_neumann = PETSC_TRUE; 163 } else { 164 is->pure_neumann = PETSC_FALSE; 165 ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 166 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 167 ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr); 168 for (i=0; i<n; i++) { 169 ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 170 } 171 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 174 ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 175 } 176 177 PetscFunctionReturn(0); 178 } 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "MatAssemblyBegin_IS" 182 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 183 { 184 Mat_IS *is = (Mat_IS*)A->data; 185 PetscErrorCode ierr; 186 187 PetscFunctionBegin; 188 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "MatAssemblyEnd_IS" 194 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 195 { 196 Mat_IS *is = (Mat_IS*)A->data; 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 EXTERN_C_BEGIN 205 #undef __FUNCT__ 206 #define __FUNCT__ "MatISGetLocalMat_IS" 207 PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat_IS(Mat mat,Mat *local) 208 { 209 Mat_IS *is = (Mat_IS *)mat->data; 210 211 PetscFunctionBegin; 212 *local = is->A; 213 PetscFunctionReturn(0); 214 } 215 EXTERN_C_END 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "MatISGetLocalMat" 219 /*@ 220 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 221 222 Input Parameter: 223 . mat - the matrix 224 225 Output Parameter: 226 . local - the local matrix usually MATSEQAIJ 227 228 Level: advanced 229 230 Notes: 231 This can be called if you have precomputed the nonzero structure of the 232 matrix and want to provide it to the inner matrix object to improve the performance 233 of the MatSetValues() operation. 234 235 .seealso: MATIS 236 @*/ 237 PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat mat,Mat *local) 238 { 239 PetscErrorCode ierr,(*f)(Mat,Mat *); 240 241 PetscFunctionBegin; 242 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 243 PetscValidPointer(local,2); 244 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr); 245 if (f) { 246 ierr = (*f)(mat,local);CHKERRQ(ierr); 247 } else { 248 local = 0; 249 } 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "MatZeroEntries_IS" 255 PetscErrorCode MatZeroEntries_IS(Mat A) 256 { 257 Mat_IS *a = (Mat_IS*)A->data; 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "MatSetOption_IS" 267 PetscErrorCode MatSetOption_IS(Mat A,MatOption op) 268 { 269 Mat_IS *a = (Mat_IS*)A->data; 270 PetscErrorCode ierr; 271 272 PetscFunctionBegin; 273 ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 277 /*MC 278 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 279 This stores the matrices in globally unassembled form. Each processor 280 assembles only its local Neumann problem and the parallel matrix vector 281 product is handled "implicitly". 282 283 Operations Provided: 284 + MatMult() 285 . MatZeroEntries() 286 . MatSetOption() 287 . MatZeroRowsLocal() 288 . MatSetValuesLocal() 289 - MatSetLocalToGlobalMapping() 290 291 Options Database Keys: 292 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 293 294 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 295 296 You must call MatSetLocalToGlobalMapping() before using this matrix type. 297 298 You can do matrix preallocation on the local matrix after you obtain it with 299 MatISGetLocalMat() 300 301 Level: advanced 302 303 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 304 305 M*/ 306 307 EXTERN_C_BEGIN 308 #undef __FUNCT__ 309 #define __FUNCT__ "MatCreate_IS" 310 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_IS(Mat A) 311 { 312 PetscErrorCode ierr; 313 Mat_IS *b; 314 315 PetscFunctionBegin; 316 ierr = PetscNew(Mat_IS,&b);CHKERRQ(ierr); 317 A->data = (void*)b; 318 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 319 A->factor = 0; 320 A->mapping = 0; 321 322 A->ops->mult = MatMult_IS; 323 A->ops->destroy = MatDestroy_IS; 324 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 325 A->ops->setvalueslocal = MatSetValuesLocal_IS; 326 A->ops->zerorowslocal = MatZeroRowsLocal_IS; 327 A->ops->assemblybegin = MatAssemblyBegin_IS; 328 A->ops->assemblyend = MatAssemblyEnd_IS; 329 A->ops->view = MatView_IS; 330 A->ops->zeroentries = MatZeroEntries_IS; 331 A->ops->setoption = MatSetOption_IS; 332 333 ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr); 334 ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr); 335 ierr = MPI_Scan(&A->m,&b->rend,1,MPIU_INT,MPI_SUM,A->comm);CHKERRQ(ierr); 336 b->rstart = b->rend - A->m; 337 338 b->A = 0; 339 b->ctx = 0; 340 b->x = 0; 341 b->y = 0; 342 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 343 344 PetscFunctionReturn(0); 345 } 346 EXTERN_C_END 347 348 349 350 351 352 353 354 355 356 357 358 359 360