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 mapping) 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,mapping,2); 169 ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); 170 if (is->mapping) { ierr = ISLocalToGlobalMappingDestroy(is->mapping);CHKERRQ(ierr); } 171 is->mapping = mapping; 172 173 /* Create the local matrix A */ 174 ierr = ISLocalToGlobalMappingGetSize(mapping,&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(mapping,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 if (rows_l) { 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 { 274 /* 275 Set up is->x as a "counting vector". This is in order to MatMult_IS 276 work properly in the interface nodes. 277 */ 278 Vec counter; 279 PetscScalar one=1.0, zero=0.0; 280 ierr = VecCreateMPI(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,&counter);CHKERRQ(ierr); 281 ierr = VecSet(counter,zero);CHKERRQ(ierr); 282 ierr = VecSet(is->x,one);CHKERRQ(ierr); 283 ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 284 ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 285 ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 286 ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 287 ierr = VecDestroy(counter);CHKERRQ(ierr); 288 } 289 if (!n) { 290 is->pure_neumann = PETSC_TRUE; 291 } else { 292 is->pure_neumann = PETSC_FALSE; 293 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 294 ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 295 for (i=0; i<n; i++) { 296 ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 297 } 298 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 299 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 301 } 302 303 PetscFunctionReturn(0); 304 } 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "MatAssemblyBegin_IS" 308 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 309 { 310 Mat_IS *is = (Mat_IS*)A->data; 311 PetscErrorCode ierr; 312 313 PetscFunctionBegin; 314 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 315 PetscFunctionReturn(0); 316 } 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "MatAssemblyEnd_IS" 320 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 321 { 322 Mat_IS *is = (Mat_IS*)A->data; 323 PetscErrorCode ierr; 324 325 PetscFunctionBegin; 326 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 330 EXTERN_C_BEGIN 331 #undef __FUNCT__ 332 #define __FUNCT__ "MatISGetLocalMat_IS" 333 PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat_IS(Mat mat,Mat *local) 334 { 335 Mat_IS *is = (Mat_IS *)mat->data; 336 337 PetscFunctionBegin; 338 *local = is->A; 339 PetscFunctionReturn(0); 340 } 341 EXTERN_C_END 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "MatISGetLocalMat" 345 /*@ 346 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 347 348 Input Parameter: 349 . mat - the matrix 350 351 Output Parameter: 352 . local - the local matrix usually MATSEQAIJ 353 354 Level: advanced 355 356 Notes: 357 This can be called if you have precomputed the nonzero structure of the 358 matrix and want to provide it to the inner matrix object to improve the performance 359 of the MatSetValues() operation. 360 361 .seealso: MATIS 362 @*/ 363 PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat mat,Mat *local) 364 { 365 PetscErrorCode ierr; 366 367 PetscFunctionBegin; 368 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 369 PetscValidPointer(local,2); 370 ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));CHKERRQ(ierr); 371 PetscFunctionReturn(0); 372 } 373 374 #undef __FUNCT__ 375 #define __FUNCT__ "MatZeroEntries_IS" 376 PetscErrorCode MatZeroEntries_IS(Mat A) 377 { 378 Mat_IS *a = (Mat_IS*)A->data; 379 PetscErrorCode ierr; 380 381 PetscFunctionBegin; 382 ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 383 PetscFunctionReturn(0); 384 } 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "MatScale_IS" 388 PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 389 { 390 Mat_IS *is = (Mat_IS*)A->data; 391 PetscErrorCode ierr; 392 393 PetscFunctionBegin; 394 ierr = MatScale(is->A,a);CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "MatGetDiagonal_IS" 400 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 401 { 402 Mat_IS *is = (Mat_IS*)A->data; 403 PetscErrorCode ierr; 404 405 PetscFunctionBegin; 406 /* get diagonal of the local matrix */ 407 ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 408 409 /* scatter diagonal back into global vector */ 410 ierr = VecSet(v,0);CHKERRQ(ierr); 411 ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 412 ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 413 PetscFunctionReturn(0); 414 } 415 416 #undef __FUNCT__ 417 #define __FUNCT__ "MatSetOption_IS" 418 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 419 { 420 Mat_IS *a = (Mat_IS*)A->data; 421 PetscErrorCode ierr; 422 423 PetscFunctionBegin; 424 ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 425 PetscFunctionReturn(0); 426 } 427 428 #undef __FUNCT__ 429 #define __FUNCT__ "MatCreateIS" 430 /*@ 431 MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 432 process but not across processes. 433 434 Input Parameters: 435 + comm - MPI communicator that will share the matrix 436 . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 437 - map - mapping that defines the global number for each local number 438 439 Output Parameter: 440 . A - the resulting matrix 441 442 Level: advanced 443 444 Notes: See MATIS for more details 445 m and n are NOT related to the size of the map, they are the size of the part of the vector owned 446 by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 447 plus the ghost points to global indices. 448 449 .seealso: MATIS, MatSetLocalToGlobalMapping() 450 @*/ 451 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 452 { 453 PetscErrorCode ierr; 454 455 PetscFunctionBegin; 456 ierr = MatCreate(comm,A);CHKERRQ(ierr); 457 ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 458 ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 459 ierr = MatSetLocalToGlobalMapping(*A,map);CHKERRQ(ierr); 460 PetscFunctionReturn(0); 461 } 462 463 /*MC 464 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 465 This stores the matrices in globally unassembled form. Each processor 466 assembles only its local Neumann problem and the parallel matrix vector 467 product is handled "implicitly". 468 469 Operations Provided: 470 + MatMult() 471 . MatMultAdd() 472 . MatMultTranspose() 473 . MatMultTransposeAdd() 474 . MatZeroEntries() 475 . MatSetOption() 476 . MatZeroRows() 477 . MatZeroRowsLocal() 478 . MatSetValues() 479 . MatSetValuesLocal() 480 . MatScale() 481 . MatGetDiagonal() 482 - MatSetLocalToGlobalMapping() 483 484 Options Database Keys: 485 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 486 487 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 488 489 You must call MatSetLocalToGlobalMapping() before using this matrix type. 490 491 You can do matrix preallocation on the local matrix after you obtain it with 492 MatISGetLocalMat() 493 494 Level: advanced 495 496 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 497 498 M*/ 499 500 EXTERN_C_BEGIN 501 #undef __FUNCT__ 502 #define __FUNCT__ "MatCreate_IS" 503 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_IS(Mat A) 504 { 505 PetscErrorCode ierr; 506 Mat_IS *b; 507 508 PetscFunctionBegin; 509 ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 510 A->data = (void*)b; 511 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 512 A->mapping = 0; 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