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