[petsc-users] Recommendations of possible preconditioner/ksp solver/snes solver for three phase flow reservoir simulator

Chukwudi Chukwudozie chdozie at gmail.com
Sat Dec 15 19:54:17 CST 2018


Hi,

I am developing a fully implicit control volume finite element black oil
reservoir simulator using petscdmplex and snes. So far, for 2 phase flow
problems (3 by 3 block for pressure, saturation and bottomhole pressure),
my simulator runs fine for small time step sizes in the order of less than
5 days. However, increasing the time step size, it doesn't converge with -3
(-5 sometimes) diverge reason. Considering that this is a fully implicit
but nonlinear problem, I expect unconditional stability for reasonable time
step sizes.

 The background DM for the 3 block problem is a composite one with the
first two being DMPLES's (bag->plexScalNode for pressure and saturation)
and the last one, bag->WellRedun, a redundant DM for the wellbore
pressures. See below.

    ierr = DMRedundantCreate(PETSC_COMM_WORLD,0
,bag->numWells,&bag->WellRedun);CHKERRQ(ierr);

    ierr =
DMCompositeCreate(PETSC_COMM_WORLD,&bag->MultiPhasePacker);CHKERRQ(ierr);

        ierr =
DMCompositeAddDM(bag->MultiPhasePacker,bag->plexScalNode);CHKERRQ(ierr);

        ierr =
DMCompositeAddDM(bag->MultiPhasePacker,bag->plexScalNode);CHKERRQ(ierr);
        ierr =
DMCompositeAddDM(bag->MultiPhasePacker,bag->WellRedun);CHKERRQ(ierr);

I have used nested matrices as shown below. In addition, I used the
fieldsplit (composite multiplicative) preconditioner with FGRES and least
square for newton solver. See below.

ierr = MatCreateNest(PETSC_COMM_WORLD,Msize,is,Msize,is,&bK[0
],K);CHKERRQ(ierr);



    ierr = SNESGetKSP(bag->snesP,&kspP);CHKERRQ(ierr);

    ierr = KSPSetTolerances(kspP,1.e-8,1.e-8
,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);

    ierr = KSPSetType(kspP,KSPFGMRES);CHKERRQ(ierr);

    ierr = KSPSetFromOptions(kspP);CHKERRQ(ierr);

    ierr = KSPGetPC(kspP,&pcP);CHKERRQ(ierr);

    ierr = PCSetType(pcP, PCFIELDSPLIT);CHKERRQ(ierr);

    ierr =
DMCompositeGetGlobalISs(bag->MultiPhasePacker,&is);CHKERRQ(ierr);

    for(i = 0; i < Msize; i++)  ierr = PCFieldSplitSetIS(pcP,NULL
,is[i]);CHKERRQ(ierr);

    ierr =
PCFieldSplitSetType(pcP,PC_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);


Any ideas on possible precoditioner options/solver options/combinations to
improve stability will be appreciated.


Chuks
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20181215/16222312/attachment.html>


More information about the petsc-users mailing list