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