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) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 212 #endif 213 ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 214 ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 215 ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 216 PetscFunctionReturn(0); 217 } 218 219 #undef ISG2LMapSetUp 220 #undef ISG2LMapApply 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "MatSetValuesLocal_IS" 224 PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 225 { 226 PetscErrorCode ierr; 227 Mat_IS *is = (Mat_IS*)A->data; 228 229 PetscFunctionBegin; 230 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "MatZeroRows_IS" 236 PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 237 { 238 Mat_IS *is = (Mat_IS*)A->data; 239 PetscInt n_l=0, *rows_l = PETSC_NULL; 240 PetscErrorCode ierr; 241 242 PetscFunctionBegin; 243 if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 244 if (n) { 245 ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr); 246 ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 247 } 248 ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 249 ierr = PetscFree(rows_l);CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "MatZeroRowsLocal_IS" 255 PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 256 { 257 Mat_IS *is = (Mat_IS*)A->data; 258 PetscErrorCode ierr; 259 PetscInt i; 260 PetscScalar *array; 261 262 PetscFunctionBegin; 263 if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 264 { 265 /* 266 Set up is->x as a "counting vector". This is in order to MatMult_IS 267 work properly in the interface nodes. 268 */ 269 Vec counter; 270 PetscScalar one=1.0, zero=0.0; 271 ierr = MatGetVecs(A,&counter,PETSC_NULL);CHKERRQ(ierr); 272 ierr = VecSet(counter,zero);CHKERRQ(ierr); 273 ierr = VecSet(is->x,one);CHKERRQ(ierr); 274 ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 275 ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 276 ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 277 ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 278 ierr = VecDestroy(&counter);CHKERRQ(ierr); 279 } 280 if (!n) { 281 is->pure_neumann = PETSC_TRUE; 282 } else { 283 is->pure_neumann = PETSC_FALSE; 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 294 PetscFunctionReturn(0); 295 } 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "MatAssemblyBegin_IS" 299 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 300 { 301 Mat_IS *is = (Mat_IS*)A->data; 302 PetscErrorCode ierr; 303 304 PetscFunctionBegin; 305 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "MatAssemblyEnd_IS" 311 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 312 { 313 Mat_IS *is = (Mat_IS*)A->data; 314 PetscErrorCode ierr; 315 316 PetscFunctionBegin; 317 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 318 PetscFunctionReturn(0); 319 } 320 321 EXTERN_C_BEGIN 322 #undef __FUNCT__ 323 #define __FUNCT__ "MatISGetLocalMat_IS" 324 PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 325 { 326 Mat_IS *is = (Mat_IS *)mat->data; 327 328 PetscFunctionBegin; 329 *local = is->A; 330 PetscFunctionReturn(0); 331 } 332 EXTERN_C_END 333 334 #undef __FUNCT__ 335 #define __FUNCT__ "MatISGetLocalMat" 336 /*@ 337 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 338 339 Input Parameter: 340 . mat - the matrix 341 342 Output Parameter: 343 . local - the local matrix usually MATSEQAIJ 344 345 Level: advanced 346 347 Notes: 348 This can be called if you have precomputed the nonzero structure of the 349 matrix and want to provide it to the inner matrix object to improve the performance 350 of the MatSetValues() operation. 351 352 .seealso: MATIS 353 @*/ 354 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 355 { 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 360 PetscValidPointer(local,2); 361 ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));CHKERRQ(ierr); 362 PetscFunctionReturn(0); 363 } 364 365 EXTERN_C_BEGIN 366 #undef __FUNCT__ 367 #define __FUNCT__ "MatISSetLocalMat_IS" 368 PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 369 { 370 Mat_IS *is = (Mat_IS *)mat->data; 371 PetscInt nrows,ncols,orows,ocols; 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 if (is->A) { 376 ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 377 ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 378 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); 379 } 380 ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 381 ierr = MatDestroy(&is->A);CHKERRQ(ierr); 382 is->A = local; 383 384 PetscFunctionReturn(0); 385 } 386 EXTERN_C_END 387 388 #undef __FUNCT__ 389 #define __FUNCT__ "MatISSetLocalMat" 390 /*@ 391 MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 392 393 Input Parameter: 394 . mat - the matrix 395 . local - the local matrix usually MATSEQAIJ 396 397 Output Parameter: 398 399 Level: advanced 400 401 Notes: 402 This can be called if you have precomputed the local matrix and 403 want to provide it to the matrix object MATIS. 404 405 .seealso: MATIS 406 @*/ 407 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 408 { 409 PetscErrorCode ierr; 410 411 PetscFunctionBegin; 412 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 413 ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 414 PetscFunctionReturn(0); 415 } 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "MatZeroEntries_IS" 419 PetscErrorCode MatZeroEntries_IS(Mat A) 420 { 421 Mat_IS *a = (Mat_IS*)A->data; 422 PetscErrorCode ierr; 423 424 PetscFunctionBegin; 425 ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNCT__ 430 #define __FUNCT__ "MatScale_IS" 431 PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 432 { 433 Mat_IS *is = (Mat_IS*)A->data; 434 PetscErrorCode ierr; 435 436 PetscFunctionBegin; 437 ierr = MatScale(is->A,a);CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 #undef __FUNCT__ 442 #define __FUNCT__ "MatGetDiagonal_IS" 443 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 444 { 445 Mat_IS *is = (Mat_IS*)A->data; 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 /* get diagonal of the local matrix */ 450 ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 451 452 /* scatter diagonal back into global vector */ 453 ierr = VecSet(v,0);CHKERRQ(ierr); 454 ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 455 ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 456 PetscFunctionReturn(0); 457 } 458 459 #undef __FUNCT__ 460 #define __FUNCT__ "MatSetOption_IS" 461 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 462 { 463 Mat_IS *a = (Mat_IS*)A->data; 464 PetscErrorCode ierr; 465 466 PetscFunctionBegin; 467 ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "MatCreateIS" 473 /*@ 474 MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 475 process but not across processes. 476 477 Input Parameters: 478 + comm - MPI communicator that will share the matrix 479 . bs - local and global block size of the matrix 480 . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 481 - map - mapping that defines the global number for each local number 482 483 Output Parameter: 484 . A - the resulting matrix 485 486 Level: advanced 487 488 Notes: See MATIS for more details 489 m and n are NOT related to the size of the map, they are the size of the part of the vector owned 490 by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 491 plus the ghost points to global indices. 492 493 .seealso: MATIS, MatSetLocalToGlobalMapping() 494 @*/ 495 PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 496 { 497 PetscErrorCode ierr; 498 499 PetscFunctionBegin; 500 ierr = MatCreate(comm,A);CHKERRQ(ierr); 501 ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 502 ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 503 ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 504 ierr = MatSetUp(*A);CHKERRQ(ierr); 505 ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 506 PetscFunctionReturn(0); 507 } 508 509 /*MC 510 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 511 This stores the matrices in globally unassembled form. Each processor 512 assembles only its local Neumann problem and the parallel matrix vector 513 product is handled "implicitly". 514 515 Operations Provided: 516 + MatMult() 517 . MatMultAdd() 518 . MatMultTranspose() 519 . MatMultTransposeAdd() 520 . MatZeroEntries() 521 . MatSetOption() 522 . MatZeroRows() 523 . MatZeroRowsLocal() 524 . MatSetValues() 525 . MatSetValuesLocal() 526 . MatScale() 527 . MatGetDiagonal() 528 - MatSetLocalToGlobalMapping() 529 530 Options Database Keys: 531 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 532 533 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 534 535 You must call MatSetLocalToGlobalMapping() before using this matrix type. 536 537 You can do matrix preallocation on the local matrix after you obtain it with 538 MatISGetLocalMat() 539 540 Level: advanced 541 542 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 543 544 M*/ 545 546 EXTERN_C_BEGIN 547 #undef __FUNCT__ 548 #define __FUNCT__ "MatCreate_IS" 549 PetscErrorCode MatCreate_IS(Mat A) 550 { 551 PetscErrorCode ierr; 552 Mat_IS *b; 553 554 PetscFunctionBegin; 555 ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 556 A->data = (void*)b; 557 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 558 559 A->ops->mult = MatMult_IS; 560 A->ops->multadd = MatMultAdd_IS; 561 A->ops->multtranspose = MatMultTranspose_IS; 562 A->ops->multtransposeadd = MatMultTransposeAdd_IS; 563 A->ops->destroy = MatDestroy_IS; 564 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 565 A->ops->setvalues = MatSetValues_IS; 566 A->ops->setvalueslocal = MatSetValuesLocal_IS; 567 A->ops->zerorows = MatZeroRows_IS; 568 A->ops->zerorowslocal = MatZeroRowsLocal_IS; 569 A->ops->assemblybegin = MatAssemblyBegin_IS; 570 A->ops->assemblyend = MatAssemblyEnd_IS; 571 A->ops->view = MatView_IS; 572 A->ops->zeroentries = MatZeroEntries_IS; 573 A->ops->scale = MatScale_IS; 574 A->ops->getdiagonal = MatGetDiagonal_IS; 575 A->ops->setoption = MatSetOption_IS; 576 577 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 578 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 579 580 b->A = 0; 581 b->ctx = 0; 582 b->x = 0; 583 b->y = 0; 584 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 585 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr); 586 ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 587 588 PetscFunctionReturn(0); 589 } 590 EXTERN_C_END 591