1 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 PetscErrorCode MatDestroy_IS(Mat A) 20 { 21 PetscErrorCode ierr; 22 Mat_IS *b = (Mat_IS*)A->data; 23 24 PetscFunctionBegin; 25 ierr = MatDestroy(&b->A);CHKERRQ(ierr); 26 ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 27 ierr = VecDestroy(&b->x);CHKERRQ(ierr); 28 ierr = VecDestroy(&b->y);CHKERRQ(ierr); 29 ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 30 ierr = PetscFree(A->data);CHKERRQ(ierr); 31 ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 32 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "MatMult_IS" 38 PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 39 { 40 PetscErrorCode ierr; 41 Mat_IS *is = (Mat_IS*)A->data; 42 PetscScalar zero = 0.0; 43 44 PetscFunctionBegin; 45 /* scatter the global vector x into the local work vector */ 46 ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 47 ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48 49 /* multiply the local matrix */ 50 ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 51 52 /* scatter product back into global memory */ 53 ierr = VecSet(y,zero);CHKERRQ(ierr); 54 ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 55 ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "MatMultAdd_IS" 61 static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 62 { 63 Vec temp_vec; 64 PetscErrorCode ierr; 65 66 PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 67 if (v3 != v2) { 68 ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 69 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 70 } else { 71 ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 72 ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 73 ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 74 ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 75 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 76 } 77 PetscFunctionReturn(0); 78 } 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "MatMultTranspose_IS" 82 PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 83 { 84 Mat_IS *is = (Mat_IS*)A->data; 85 PetscErrorCode ierr; 86 87 PetscFunctionBegin; /* y = A' * x */ 88 /* scatter the global vector x into the local work vector */ 89 ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 90 ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 91 92 /* multiply the local matrix */ 93 ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 94 95 /* scatter product back into global vector */ 96 ierr = VecSet(y,0);CHKERRQ(ierr); 97 ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 98 ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "MatMultTransposeAdd_IS" 104 PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 105 { 106 Vec temp_vec; 107 PetscErrorCode ierr; 108 109 PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 110 if (v3 != v2) { 111 ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 112 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 113 } else { 114 ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 115 ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 116 ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 117 ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 118 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 119 } 120 PetscFunctionReturn(0); 121 } 122 123 #undef __FUNCT__ 124 #define __FUNCT__ "MatView_IS" 125 PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 126 { 127 Mat_IS *a = (Mat_IS*)A->data; 128 PetscErrorCode ierr; 129 PetscViewer sviewer; 130 131 PetscFunctionBegin; 132 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 133 ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 134 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 135 PetscFunctionReturn(0); 136 } 137 138 #undef __FUNCT__ 139 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 140 PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 141 { 142 PetscErrorCode ierr; 143 PetscInt n,bs; 144 Mat_IS *is = (Mat_IS*)A->data; 145 IS from,to; 146 Vec global; 147 148 PetscFunctionBegin; 149 if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix"); 150 PetscCheckSameComm(A,1,rmapping,2); 151 if (rmapping != cmapping) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 152 ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 153 ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 154 is->mapping = rmapping; 155 156 /* Create the local matrix A */ 157 ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 158 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 159 ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 160 ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 161 ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 162 ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 163 ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 164 ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 165 166 /* Create the local work vectors */ 167 ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 168 ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 169 ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 170 ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 171 ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 172 ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 173 ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 174 175 /* setup the global to local scatter */ 176 ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 177 ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 178 ierr = MatGetVecs(A,&global,PETSC_NULL);CHKERRQ(ierr); 179 ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 180 ierr = VecDestroy(&global);CHKERRQ(ierr); 181 ierr = ISDestroy(&to);CHKERRQ(ierr); 182 ierr = ISDestroy(&from);CHKERRQ(ierr); 183 PetscFunctionReturn(0); 184 } 185 186 #define ISG2LMapApply(mapping,n,in,out) 0; \ 187 if (!(mapping)->globals) { \ 188 PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr); \ 189 } \ 190 { \ 191 PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend; \ 192 for (_i=0; _i<n; _i++) { \ 193 if (in[_i] < 0) out[_i] = in[_i]; \ 194 else if (in[_i] < _start) out[_i] = -1; \ 195 else if (in[_i] > _end) out[_i] = -1; \ 196 else out[_i] = _globals[in[_i] - _start]; \ 197 } \ 198 } 199 200 #undef __FUNCT__ 201 #define __FUNCT__ "MatSetValues_IS" 202 PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 203 { 204 Mat_IS *is = (Mat_IS*)mat->data; 205 PetscInt rows_l[2048],cols_l[2048]; 206 PetscErrorCode ierr; 207 208 PetscFunctionBegin; 209 #if defined(PETSC_USE_DEBUG) 210 if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 211 #endif 212 ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 213 ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 214 ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 215 PetscFunctionReturn(0); 216 } 217 218 #undef ISG2LMapSetUp 219 #undef ISG2LMapApply 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "MatSetValuesLocal_IS" 223 PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 224 { 225 PetscErrorCode ierr; 226 Mat_IS *is = (Mat_IS*)A->data; 227 228 PetscFunctionBegin; 229 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232 233 #undef __FUNCT__ 234 #define __FUNCT__ "MatZeroRows_IS" 235 PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 236 { 237 Mat_IS *is = (Mat_IS*)A->data; 238 PetscInt n_l = 0, *rows_l = PETSC_NULL; 239 PetscErrorCode ierr; 240 241 PetscFunctionBegin; 242 if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 243 if (n) { 244 ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr); 245 ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 246 } 247 ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 248 ierr = PetscFree(rows_l);CHKERRQ(ierr); 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "MatZeroRowsLocal_IS" 254 PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 255 { 256 Mat_IS *is = (Mat_IS*)A->data; 257 PetscErrorCode ierr; 258 PetscInt i; 259 PetscScalar *array; 260 261 PetscFunctionBegin; 262 if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 263 { 264 /* 265 Set up is->x as a "counting vector". This is in order to MatMult_IS 266 work properly in the interface nodes. 267 */ 268 Vec counter; 269 PetscScalar one=1.0, zero=0.0; 270 ierr = MatGetVecs(A,&counter,PETSC_NULL);CHKERRQ(ierr); 271 ierr = VecSet(counter,zero);CHKERRQ(ierr); 272 ierr = VecSet(is->x,one);CHKERRQ(ierr); 273 ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 274 ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 275 ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 276 ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 277 ierr = VecDestroy(&counter);CHKERRQ(ierr); 278 } 279 if (!n) { 280 is->pure_neumann = PETSC_TRUE; 281 } else { 282 is->pure_neumann = PETSC_FALSE; 283 284 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 285 ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 286 for (i=0; i<n; i++) { 287 ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 288 } 289 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 290 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 291 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 292 } 293 PetscFunctionReturn(0); 294 } 295 296 #undef __FUNCT__ 297 #define __FUNCT__ "MatAssemblyBegin_IS" 298 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 299 { 300 Mat_IS *is = (Mat_IS*)A->data; 301 PetscErrorCode ierr; 302 303 PetscFunctionBegin; 304 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 305 PetscFunctionReturn(0); 306 } 307 308 #undef __FUNCT__ 309 #define __FUNCT__ "MatAssemblyEnd_IS" 310 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 311 { 312 Mat_IS *is = (Mat_IS*)A->data; 313 PetscErrorCode ierr; 314 315 PetscFunctionBegin; 316 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 320 EXTERN_C_BEGIN 321 #undef __FUNCT__ 322 #define __FUNCT__ "MatISGetLocalMat_IS" 323 PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 324 { 325 Mat_IS *is = (Mat_IS*)mat->data; 326 327 PetscFunctionBegin; 328 *local = is->A; 329 PetscFunctionReturn(0); 330 } 331 EXTERN_C_END 332 333 #undef __FUNCT__ 334 #define __FUNCT__ "MatISGetLocalMat" 335 /*@ 336 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 337 338 Input Parameter: 339 . mat - the matrix 340 341 Output Parameter: 342 . local - the local matrix usually MATSEQAIJ 343 344 Level: advanced 345 346 Notes: 347 This can be called if you have precomputed the nonzero structure of the 348 matrix and want to provide it to the inner matrix object to improve the performance 349 of the MatSetValues() operation. 350 351 .seealso: MATIS 352 @*/ 353 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 354 { 355 PetscErrorCode ierr; 356 357 PetscFunctionBegin; 358 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 359 PetscValidPointer(local,2); 360 ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 EXTERN_C_BEGIN 365 #undef __FUNCT__ 366 #define __FUNCT__ "MatISSetLocalMat_IS" 367 PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 368 { 369 Mat_IS *is = (Mat_IS*)mat->data; 370 PetscInt nrows,ncols,orows,ocols; 371 PetscErrorCode ierr; 372 373 PetscFunctionBegin; 374 if (is->A) { 375 ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 376 ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 377 if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols); 378 } 379 ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 380 ierr = MatDestroy(&is->A);CHKERRQ(ierr); 381 is->A = local; 382 PetscFunctionReturn(0); 383 } 384 EXTERN_C_END 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "MatISSetLocalMat" 388 /*@ 389 MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 390 391 Input Parameter: 392 . mat - the matrix 393 . local - the local matrix usually MATSEQAIJ 394 395 Output Parameter: 396 397 Level: advanced 398 399 Notes: 400 This can be called if you have precomputed the local matrix and 401 want to provide it to the matrix object MATIS. 402 403 .seealso: MATIS 404 @*/ 405 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 406 { 407 PetscErrorCode ierr; 408 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 411 ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 412 PetscFunctionReturn(0); 413 } 414 415 #undef __FUNCT__ 416 #define __FUNCT__ "MatZeroEntries_IS" 417 PetscErrorCode MatZeroEntries_IS(Mat A) 418 { 419 Mat_IS *a = (Mat_IS*)A->data; 420 PetscErrorCode ierr; 421 422 PetscFunctionBegin; 423 ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 424 PetscFunctionReturn(0); 425 } 426 427 #undef __FUNCT__ 428 #define __FUNCT__ "MatScale_IS" 429 PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 430 { 431 Mat_IS *is = (Mat_IS*)A->data; 432 PetscErrorCode ierr; 433 434 PetscFunctionBegin; 435 ierr = MatScale(is->A,a);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "MatGetDiagonal_IS" 441 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 442 { 443 Mat_IS *is = (Mat_IS*)A->data; 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 /* get diagonal of the local matrix */ 448 ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 449 450 /* scatter diagonal back into global vector */ 451 ierr = VecSet(v,0);CHKERRQ(ierr); 452 ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 453 ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 454 PetscFunctionReturn(0); 455 } 456 457 #undef __FUNCT__ 458 #define __FUNCT__ "MatSetOption_IS" 459 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 460 { 461 Mat_IS *a = (Mat_IS*)A->data; 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "MatCreateIS" 471 /*@ 472 MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 473 process but not across processes. 474 475 Input Parameters: 476 + comm - MPI communicator that will share the matrix 477 . bs - local and global block size of the matrix 478 . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 479 - map - mapping that defines the global number for each local number 480 481 Output Parameter: 482 . A - the resulting matrix 483 484 Level: advanced 485 486 Notes: See MATIS for more details 487 m and n are NOT related to the size of the map, they are the size of the part of the vector owned 488 by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 489 plus the ghost points to global indices. 490 491 .seealso: MATIS, MatSetLocalToGlobalMapping() 492 @*/ 493 PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 494 { 495 PetscErrorCode ierr; 496 497 PetscFunctionBegin; 498 ierr = MatCreate(comm,A);CHKERRQ(ierr); 499 ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 500 ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 501 ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 502 ierr = MatSetUp(*A);CHKERRQ(ierr); 503 ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 504 PetscFunctionReturn(0); 505 } 506 507 /*MC 508 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 509 This stores the matrices in globally unassembled form. Each processor 510 assembles only its local Neumann problem and the parallel matrix vector 511 product is handled "implicitly". 512 513 Operations Provided: 514 + MatMult() 515 . MatMultAdd() 516 . MatMultTranspose() 517 . MatMultTransposeAdd() 518 . MatZeroEntries() 519 . MatSetOption() 520 . MatZeroRows() 521 . MatZeroRowsLocal() 522 . MatSetValues() 523 . MatSetValuesLocal() 524 . MatScale() 525 . MatGetDiagonal() 526 - MatSetLocalToGlobalMapping() 527 528 Options Database Keys: 529 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 530 531 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 532 533 You must call MatSetLocalToGlobalMapping() before using this matrix type. 534 535 You can do matrix preallocation on the local matrix after you obtain it with 536 MatISGetLocalMat() 537 538 Level: advanced 539 540 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 541 542 M*/ 543 544 EXTERN_C_BEGIN 545 #undef __FUNCT__ 546 #define __FUNCT__ "MatCreate_IS" 547 PetscErrorCode MatCreate_IS(Mat A) 548 { 549 PetscErrorCode ierr; 550 Mat_IS *b; 551 552 PetscFunctionBegin; 553 ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 554 A->data = (void*)b; 555 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 556 557 A->ops->mult = MatMult_IS; 558 A->ops->multadd = MatMultAdd_IS; 559 A->ops->multtranspose = MatMultTranspose_IS; 560 A->ops->multtransposeadd = MatMultTransposeAdd_IS; 561 A->ops->destroy = MatDestroy_IS; 562 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 563 A->ops->setvalues = MatSetValues_IS; 564 A->ops->setvalueslocal = MatSetValuesLocal_IS; 565 A->ops->zerorows = MatZeroRows_IS; 566 A->ops->zerorowslocal = MatZeroRowsLocal_IS; 567 A->ops->assemblybegin = MatAssemblyBegin_IS; 568 A->ops->assemblyend = MatAssemblyEnd_IS; 569 A->ops->view = MatView_IS; 570 A->ops->zeroentries = MatZeroEntries_IS; 571 A->ops->scale = MatScale_IS; 572 A->ops->getdiagonal = MatGetDiagonal_IS; 573 A->ops->setoption = MatSetOption_IS; 574 575 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 576 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 577 578 b->A = 0; 579 b->ctx = 0; 580 b->x = 0; 581 b->y = 0; 582 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 583 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr); 584 ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 585 PetscFunctionReturn(0); 586 } 587 EXTERN_C_END 588