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