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