[petsc-dev] KSPsolve

Mark Adams mfadams at lbl.gov
Sat Oct 3 04:51:55 CDT 2015


Noel,

I am cc'ing the PETSc mailing list. BTW, this question is very normal and,
in hindsight, you could have sent us this earlier.

Plasticity can be very hard, but you can do a lot better than ILU.  Your
examples so far do not look terrible, but plasticity can get really hard
(or stiff to be more precise).  First, your couple thermal problems will
probably want to use FieldSplit, but let's work on solid mechanics part
first.

Are these symmetric positive definite (SPD, normal finite element
discretizations)?  I assume they are, or at least are close enough.  If
they are then you want to use -ksp_type cg. I've appended a decent
parameter list to start; this will get tweaked with optimizations,
depending on your problems.  You need to use AMG and you need to give AMG
the null space of the operator.  See
petsc/src/ksp/ksp/example/tutorial/ex56.c for an example of a 3D elasticity
problem.  You need to give us coordinates of the vertices basically.

Mark

-ksp_type cg
-ksp_rtol 1.e-6
-pc_type gamg
-pc_gamg_agg_nsmooths 1
-pc_gamg_threshold -0.02
-gamg_est_ksp_type cg
#-mg_levels_ksp_type richardson
-mg_levels_ksp_type chebyshev
-mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05
-mg_levels_pc_type sor
-mg_levels_ksp_max_it 2
-ksp_converged_reason
#-ksp_view
-pc_gamg_repartition true
-pc_gamg_mat_partitioning_type parmetis
-pc_gamg_process_eq_limit 200
#-pc_mg_log
-ksp_monitor


On Fri, Oct 2, 2015 at 7:02 PM, Noel Keen <ndkeen at lbl.gov> wrote:

> Hi Mark,
>
> I don't mind admitting that I'm not very familiar with the physical
> problem and the equations,
> so I had to ask someone else.  I was hoping it would be a simple matter of
> replacing the serial
> Ax=b with a parallel one, but it's probably never that easy for the
> interesting problems.  :)
>
> From others (T+M stands for Tough+Mechanics or the Tough fluid code
> coupled with a geomechanics code [named rocmech])
> 1. We have both 2D and 3D models for geomechanics (solid model).
> 2. Plasticity is there. We have several publications for plasticity of
> T+M. But, regarding the hydrate case, we have not tried yet, although it
> certainly works.
> (I can verify that I see this in the code and we do *not* use it for these
> inputs I'm trying to work with)
> 3. We have poroelasticity. Precisely, we have rigorously two-way coupled
> thermo-poro-elasto-plasticity. T+M is the most advanced simulator in
> coupled flow and geomechanics. People in this community already know this,
> I believe.
> 4. The most incompressible fluid is water (almost cw=4.4e-10 Pa^{-1}).
> Although civil engineering people approximate water compressibility to
> zero,  the accurate one is TOUGH. Also, if we deal with multiphase flow,
> flow is certainly compressible. e.g., for oil-water phase, ct(total fluid
> compressibility)=Soco+Swcw, sum of each saturation times compressibility.
> (That is, depending on the EOS, the fluid could be compressible)
>
>
> Noel
>
> On Oct 2, 2015, at 12:52 PM, Mark Adams wrote:
>
> > Noel,
> >
> > It sounds like you are doing 3D solids.  Geomechanics problems use
> plasticity a lot, which is hard.  I don't know what you doing.
> >
> > Are these Stokes problems or are they compressible?
> >
> > Mark
> >
> >
> > On Fri, Oct 2, 2015 at 12:32 PM, noel keen LBL <ndkeen at lbl.gov> wrote:
> >
> >
> > Begin forwarded message:
> >
> > > From: noel keen LBL <ndkeen at lbl.gov>
> > > Date: October 2, 2015 10:31:22 AM PDT
> > > To: "Mark F. Adams" <mark.adams at columbia.edu>, noel keen <
> ndkeen at lbl.gov>
> > > Subject: KSPsolve
> > >
> > > Hello Mark,
> > >
> > > At on point I asked if you could help with PETSc.  I was able to
> figure it out myself following examples, etc.
> > > And then sidetracked with other tasks, etc.
> > >
> > > I've finally reached a point where the bottleneck is the PETSc
> parallel solve in our code (tough+rocmech http://tinyurl.com/nskrgf8).
> > > There are essentially 2 "types" of solves, one is the standard tough
> "flow" and the other is the newer geomechanical solve (aka the rocmech
> solve).
> > > I'm having trouble with both, but certainly the rocmech solve is the
> slowest.  I created a stand-alone PETSc code (written in F90,
> > > unfortunately), that simply reads in a Jacobian (which is output from
> an actual simulation run) and then solves the system of equations
> > > using KSPsolve.  I can then solve the system again using their serial
> solver and compare the performance and results directly.
> > >
> > > The serial solver has been surprisingly robust and efficient for them
> http://tinyurl.com/opgjgw5.  It uses a version of ILU and CGS.
> > > I was convinced (and still am), that once I got PETSc setup properly,
> I would be able to obtain a large speedup.  I can just beat it in
> > > serial, but when I add procs, PETSc really struggles.  I've tried
> hundreds if not thousands of different permutations of PC's, PC options,
> > > KSP options, and numbers of processors.  I _can_ get a speedup, but am
> just baffled at why the performance (and stability) is not better.
> > >
> > > Can you offer any help?  I can show/give you my small stand-alone code
> and the data.
> > > Could you recommend someone else to help?
> > > Should I even expect a parallel solver to significantly outperform an
> efficient serial solve?
> > > Are there any obvious "tricks" the serial solve is doing that would
> help us in determining the best way to choose parallel solver?
> > >
> > > Some basic stats and sample of results:
> > >
> > > nrows=      261144  nnzero=    19879452
> > >
> > > I have the preallocation done correctly (I know this will destroy
> parallel perf) . For a parallel run, all procs report:
> > > MatGetInfo says nnz=      19879452. nlocal=     21762 nrows=    261144
> Mallocs=          0. NZalloc=    19879452. NZunneeded=           0.
> > >
> > > For this system, the serial solver can reach a tolerance of 1.0e-6 in
> about 22 seconds on edison.  So that's the time to beat.
> > >
> > > Some of the problems I'm seeing:  poor performance, method divergence
> as I decrease tolerance or increase number of procs,
> > > num iterations/performance don't always make sense when changing
> options/parameters, number of iterations often
> > > *increasing* as I add procs (always bothersome).
> > >
> > > I typically test with a KSP tolerance of 1.0e-4 as I try options.
> Here are some good results:
> > >
> > > -ksp_type bcgs -ksp_rtol 1.e-4 -pc_type bjacobi -sub_pc_type ilu
> > > kspsolve np=  48 num_iterations=     562 rnorm=          4.0257E-10
> reason=   2 wcsolve=    2.56
> > > kspsolve np=  24 num_iterations=     470 rnorm=          6.3274E-10
> reason=   2 wcsolve=    2.42
> > > kspsolve np=  12 num_iterations=     238 rnorm=          7.4453E-10
> reason=   2 wcsolve=    2.08
> > >
> > > -ksp_view ascii -ksp_type bcgs -ksp_rtol 1.e-4 -pc_type hypre
> -pc_hypre_type pilut -pc_hypre_pilut_maxiter 200 -pc_hypre_pilut_tol 1.0e-4
> > > kspsolve np=  96 num_iterations=    1053 rnorm=          2.7155E-10
> reason=   2 wcsolve=   15.16
> > > kspsolve np=  48 num_iterations=     913 rnorm=          2.8016E-10
> reason=   2 wcsolve=    1.64
> > > kspsolve np=  24 num_iterations=     794 rnorm=          2.7149E-10
> reason=   2 wcsolve=    2.75
> > > kspsolve np=   4 num_iterations=    2375 rnorm=          3.0389E-10
> reason=   2 wcsolve=   38.19
> > >
> > >
> > > But then going to 1.0e-6 KSP tolerance (a fair comparison to serial
> DLUST solver) I have trouble:
> > >
> > > -ksp_view ascii -ksp_type bcgs -ksp_rtol 1.e-6 -pc_type bjacobi
> -sub_pc_type ilu
> > > kspsolve np=  48 num_iterations=    2179 rnorm=          3.0013E-12
> reason=   2 wcsolve=    9.26
> > > kspsolve np=  24 num_iterations=    1100 rnorm=          2.2589E-10
> reason=  -5 wcsolve=    5.63  <-- diverges
> > > kspsolve np=  12 num_iterations=     759 rnorm=          4.9629E-07
> reason=  -5 wcsolve=    6.38 <-- diverges
> > >
> > > -ksp_view ascii -ksp_type bcgs -ksp_rtol 1.e-6 -pc_type hypre
> -pc_hypre_type pilut -pc_hypre_pilut_maxitr 1000 -pc_hypre_pilut_tol 1.0e-10
> > > kspsolve np=  48 num_iterations=    7579 rnorm=          5.7098E-12
> reason=   2 wcsolve=   43.81 <-- slow, lots of iterations
> > >
> > >
> > >
> > > Thanks
> > >
> > > Noel
> >
> >
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20151003/55804129/attachment.html>


More information about the petsc-dev mailing list