[petsc-users] [DMMG] Stokes Solver

domenico.borzacchiello at univ-st-etienne.fr domenico.borzacchiello at univ-st-etienne.fr
Fri Apr 8 04:40:13 CDT 2011


Hi Jed
Thank you for the very quick reply,

running with -ksp_converged_reason and lu preconditioner gives out

Linear solve converged due to CONVERGED_RTOL iterations 14
Linear solve converged due to CONVERGED_RTOL iterations 12
Linear solve converged due to CONVERGED_RTOL iterations 13
Linear solve converged due to CONVERGED_RTOL iterations 12
Linear solve converged due to CONVERGED_RTOL iterations 13
Linear solve converged due to CONVERGED_RTOL iterations 12
Linear solve converged due to CONVERGED_RTOL iterations 12
Linear solve did not converge due to DIVERGED_DTOL iterations 1080
      RINFO(1) (local estimated flops for the elimination after analysis):
             [0] 5.42609e+08
      RINFO(2) (local estimated flops for the assembly after factorization):
             [0]  3.83582e+06
      RINFO(3) (local estimated flops for the elimination after
factorization):
             [0]  5.44162e+08
      INFO(15) (estimated size of (in MB) MUMPS internal data for running
numerical factorization):
             [0] 17
      INFO(16) (size of (in MB) MUMPS internal data used during numerical
factorization):
             [0] 17
      INFO(23) (num of pivots eliminated on this processor after
factorization):
             [0] 1372
Number of Newton iterations = 7
Number of Linear iterations = 88
Average Linear its / Newton = 1.257143e+01
Converged Reason = -3

I can't get why it takes more than one Newton iteration even if my system
is linear and the linear solver is direct.

As for the Vanka smoother would you suggest to implement it with PCSHELL
or by defining a new PC type?

thanks,
Domenico.

> On Fri, Apr 8, 2011 at 09:59,
> <domenico.borzacchiello at univ-st-etienne.fr>wrote:
>
>> Hi,
>> I'm trying to implement a 3d stokes solver on a simple cartesian
>> staggered
>> grid using the a MAC FV Discretisation. I'm following the example given
>> in
>> /snes/examples/tutorials/ex30.c. As the stokes solver is ready I'll add
>> some more complex constitutive laws. For the moment I'm testing the
>> solver
>> with just 1-level MG and I'd like a few clarification to know if what
>> I'm
>> doing is correct.
>>
>> - The solver solves for u v w p. I'm using a single DA with 4 DOFs and
>> due
>> to the MAC arrangement I have some "extraboundary" nodes for u ( in the
>> y
>> and z directions) v (x & z dir) w ( x & y dir) and p (x,y & z dir). If
>> what I'm getting from ex30.c is right I have to write a simple identity
>> for each of these nodes (i.e. p_extra = anyvalue) as they are not
>> coupled
>> with the rest of the system. I'm doing the same for Dirichlet BCs nodes
>> (i.e. u = Ubound). Is this correct?
>>
>
> yes
>
>
>>
>> - How does Petsc deal with the pressure-velocity coupling? Is it correct
>> to try to solve the whole coupled system with DMMG as in ex30.c? At
>> present time I'm getting no convergence by running dmmgsolve (snes) with
>> all the default options on a very small system.
>>
>>  0 SNES Function norm 7.128085632250e+00
>> Number of Newton iterations = 0
>> Number of Linear iterations = 0
>> Average Linear its / Newton = -nan
>> Converged Reason = -3
>>
>
> You may as well check why the linear solve failed by running with
> -ksp_converged_reason.
>
> There are two challenges for solving Stokes in this way. First, the
> interpolation operators that the DA gives you are probably not what you
> want
> (pressure and velocity should be interpolated differently) and second, the
> standard smoother of ILU is not expected to work with interlaced velocity
> and pressure (you would want to either use a "Vanka smoother" that solves
> small problems associated with each pressure cell and all adjacent
> velocities, use field-split as a smoother, or (tricky and not guaranteed
> to
> work) order the pressures last in each subdomain with ILU as a smoother).
>
> Vanka smoothers are very problem-dependent so you would need to write that
> part yourself. It would be nice to have an example for Stokes. An
> alternative to coupled multigrid is to use PCFieldSplit at the top level
> and
> multigrid for the viscous part separately (see e.g. Elman's many papers on
> this approach). We don't currently have an example doing PCFieldSplit with
> geometric multigrid inside the splits, but it should work with petsc-dev
> if
> you follow the approach in src/ksp/ksp/examples/tutorials/ex45.c (which
> does
> not use DMMG, we are working to fold DMMG's functionality into the
> existing
> solver objects).
>




More information about the petsc-users mailing list