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