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