Changeset 23925
- Timestamp:
- 05/20/19 16:33:29 (6 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/solutionsequences/solutionsequence_schurcg.cpp
r23924 r23925 88 88 * and I the Schur preconditioner (stored here, because the space was allocated either way) 89 89 * */ 90 #if _PETSC_MINOR_>10 91 MatCreateSubMatrix(Kff,isv,isv,MAT_INITIAL_MATRIX,&A); 92 MatCreateSubMatrix(Kff,isv,isp,MAT_INITIAL_MATRIX,&B); 93 MatCreateSubMatrix(Kff,isp,isv,MAT_INITIAL_MATRIX,&BT); 94 #else 90 95 MatGetSubMatrix(Kff,isv,isv,MAT_INITIAL_MATRIX,&A); 91 96 MatGetSubMatrix(Kff,isv,isp,MAT_INITIAL_MATRIX,&B); 92 97 MatGetSubMatrix(Kff,isp,isv,MAT_INITIAL_MATRIX,&BT); 98 #endif 93 99 94 100 /* Extract preconditioner matrix on the pressure space*/ 101 #if _PETSC_MINOR_>10 102 MatCreateSubMatrix(Kff,isp,isp,MAT_INITIAL_MATRIX,&IP); 103 #else 95 104 MatGetSubMatrix(Kff,isp,isp,MAT_INITIAL_MATRIX,&IP); 105 #endif 96 106 97 107 … … 522 532 } 523 533 524 /* Note: SchurCG also constructs the Schur preconditioner and stores it in the free block of Kff */525 /* [A B]534 /* Note: SchurCG also constructs the Schur preconditioner and stores it in the free block of Kff 535 * [A B] 526 536 * Kff = | | 527 537 * [B^T IP] 528 /* To calculate the residual, only the necessary blocks need to be extracted */ 529 530 /*Extract A, B, B^T */ 531 MatGetSubMatrix(Kff->pmatrix->matrix,isv,isv,MAT_INITIAL_MATRIX,&A); 532 MatGetSubMatrix(Kff->pmatrix->matrix,isv,isp,MAT_INITIAL_MATRIX,&B); 533 MatGetSubMatrix(Kff->pmatrix->matrix,isp,isv,MAT_INITIAL_MATRIX,&BT); 538 * To calculate the residual, only the necessary blocks need to be extracted */ 539 540 /*Extract A, B, B^T */ 541 #if _PETSC_MINOR_>10 542 MatCreateSubMatrix(Kff->pmatrix->matrix,isv,isv,MAT_INITIAL_MATRIX,&A); 543 MatCreateSubMatrix(Kff->pmatrix->matrix,isv,isp,MAT_INITIAL_MATRIX,&B); 544 MatCreateSubMatrix(Kff->pmatrix->matrix,isp,isv,MAT_INITIAL_MATRIX,&BT); 545 #else 546 MatGetSubMatrix(Kff->pmatrix->matrix,isv,isv,MAT_INITIAL_MATRIX,&A); 547 MatGetSubMatrix(Kff->pmatrix->matrix,isv,isp,MAT_INITIAL_MATRIX,&B); 548 MatGetSubMatrix(Kff->pmatrix->matrix,isp,isv,MAT_INITIAL_MATRIX,&BT); 549 #endif 534 550 535 551 /*no. of free nodes in velocity/pressure space*/ … … 688 704 Vector<IssmDouble>* df = NULL; 689 705 Vector<IssmDouble>* ys = NULL; 690 691 692 706 693 707 /*parameters:*/ … … 707 721 int count=0; 708 722 bool converged=false; 709 710 723 711 724 /*Start non-linear iteration using input velocity: */ … … 716 729 InputUpdateFromConstantx(femmodel,converged,ConvergedEnum); 717 730 InputUpdateFromSolutionx(femmodel,ug); 718 719 731 720 732 for(;;){ … … 724 736 delete ug; 725 737 726 727 728 738 /*Get stiffness matrix and Load vector*/ 729 739 SystemMatricesx(&Kff,&Kfs,&pf,&df,NULL,femmodel);
Note:
See TracChangeset
for help on using the changeset viewer.