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