Changeset 23925


Ignore:
Timestamp:
05/20/19 16:33:29 (6 years ago)
Author:
Mathieu Morlighem
Message:

CHG: fixing warnings with new version of PETSc

File:
1 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk-jpl/src/c/solutionsequences/solutionsequence_schurcg.cpp

    r23924 r23925  
    8888         * and I the Schur preconditioner (stored here, because the space was allocated either way)
    8989         *         */
     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
    9095        MatGetSubMatrix(Kff,isv,isv,MAT_INITIAL_MATRIX,&A);
    9196        MatGetSubMatrix(Kff,isv,isp,MAT_INITIAL_MATRIX,&B);
    9297        MatGetSubMatrix(Kff,isp,isv,MAT_INITIAL_MATRIX,&BT);
     98        #endif
    9399       
    94100        /* Extract preconditioner matrix on the pressure space*/
     101        #if _PETSC_MINOR_>10
     102        MatCreateSubMatrix(Kff,isp,isp,MAT_INITIAL_MATRIX,&IP);
     103        #else
    95104        MatGetSubMatrix(Kff,isp,isp,MAT_INITIAL_MATRIX,&IP);
     105        #endif
    96106
    97107
     
    522532        }
    523533
    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]
    526536        * Kff =  |      |
    527537        *                       [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
    534550       
    535551                /*no. of free nodes in velocity/pressure space*/
     
    688704        Vector<IssmDouble>* df  = NULL;
    689705        Vector<IssmDouble>* ys  = NULL;
    690 
    691 
    692706
    693707        /*parameters:*/
     
    707721        int  count=0;
    708722        bool converged=false;
    709 
    710723       
    711724        /*Start non-linear iteration using input velocity: */
     
    716729        InputUpdateFromConstantx(femmodel,converged,ConvergedEnum);
    717730        InputUpdateFromSolutionx(femmodel,ug);
    718 
    719731       
    720732        for(;;){
     
    724736                delete ug;
    725737
    726 
    727                
    728738                /*Get stiffness matrix and Load vector*/
    729739                SystemMatricesx(&Kff,&Kfs,&pf,&df,NULL,femmodel);
Note: See TracChangeset for help on using the changeset viewer.