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