HYPRE with multiple variables

Nicolas Bathfield nicolas.bathfield at chalmers.se
Tue May 8 04:32:34 CDT 2007


Hi,


Boomeramg is used as preconditioner, and I use default gmres to solve the
equations. I did not expect the preconditioner to be used at every
iteration, but I now understand better what's going on.
I solve a spase matrix with 5 diagonals.
I tried to use boomeramg as a solver only. To achieve this, I had to
comment out the command  ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE),
and I used the option -ksp_preonly. Is it the right way to do?

Best regards,

Nicolas

> Hi,
>
> You mentioned that your preconditioner is boomerAMG. What solver are
> you using?
> What type of matrix are you trying to solve? Is it SPD and sparse?
> Either way, if your solver is not AMG and you are using AMG as a
> preconditioner, I would suggest that you use a stopping criteria for
> AMG such as a relative residual convergence of say 1e-03 for the
> preconditioner, and the solvers relative convergence 1e-05 as in your
> example, and disregard a stopping criteria which involves maximum
> number of iterations for AMG, set that to something high. Maybe you
> will get better results that way. In the end, I guess it is a lot of
> testing and fine-tuning.
>
> With best regards, Shaman Mahmoudi
>
> On May 7, 2007, at 4:50 PM, Nicolas Bathfield wrote:
>
>> Hi Allison, Barry, and PETSc users,
>>
>> The nodal relaxation leads to a nice solution! Thanks a lot for
>> your help
>> on this!
>>
>> Even though I get a satisfactory result, I can not explain why
>> hypre is
>> performing what looks like several outer loops (see below for an
>> abstract
>> of the output from -pc_hypre_boomeramg_print_statistics).
>>
>> Here is the list of options I used:
>>
>> -pc_hypre_boomeramg_tol 1e-5
>> -pc_type hypre -pc_hypre_type boomeramg
>> -pc_hypre_boomeramg_print_statistics
>> -pc_hypre_boomeramg_grid_sweeps_all 20
>> -pc_hypre_boomeramg_max_iter 6
>> -pc_hypre_boomeramg_relax_weight_all 0.6
>> -pc_hypre_boomeramg_outer_relax_weight_all 0.6
>> -pc_hypre_boomeramg_nodal_coarsen 1
>>
>> I would expect the preconditioner to be executed once only (in my case
>> this includes 6 V cycles).
>> As you can see from the print_statistic option, boomeramg is exectuted
>> several times (many more times than displyed here, it was too long
>> to copy
>> everything). Is this normal? If yes, why is it so?
>>
>> Best regards,
>>
>> Nicolas
>>
>> BoomerAMG SOLVER PARAMETERS:
>>
>>   Maximum number of cycles:         6
>>   Stopping Tolerance:               1.000000e-05
>>   Cycle type (1 = V, 2 = W, etc.):  1
>>
>>   Relaxation Parameters:
>>    Visiting Grid:                     down   up  coarse
>>    Visiting Grid:                     down   up  coarse
>>             Number of partial sweeps:   20   20     1
>>    Type 0=Jac, 3=hGS, 6=hSGS, 9=GE:      6    6     6
>>    Point types, partial sweeps (1=C, -1=F):
>>                   Pre-CG relaxation (down):   1  -1   1  -1   1
>> -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1
>>                    Post-CG relaxation (up):  -1   1  -1   1  -1
>> 1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1
>>                              Coarsest grid:   0
>>
>>  Relaxation Weight 0.600000 level 0
>>  Outer relaxation weight 0.600000 level 0
>>  Output flag (print_level): 3
>>
>>
>> AMG SOLUTION INFO:
>>                                             relative
>>                residual        factor       residual
>>                --------        ------       --------
>>     Initial    4.430431e-03                 1.000000e+00
>>     Cycle  1   4.540275e-03    1.024793     1.024793e+00
>>     Cycle  2   4.539148e-03    0.999752     1.024539e+00
>>     Cycle  3   4.620010e-03    1.017814     1.042790e+00
>>     Cycle  4   5.196532e-03    1.124788     1.172918e+00
>>     Cycle  5   6.243043e-03    1.201386     1.409128e+00
>>     Cycle  6   7.310008e-03    1.170905     1.649954e+00
>>
>>
>> ==============================================
>>  NOTE: Convergence tolerance was not achieved
>>       within the allowed 6 V-cycles
>> ==============================================
>>
>>  Average Convergence Factor = 1.087039
>>
>>      Complexity:    grid = 1.000000
>>                 operator = 1.000000
>>                    cycle = 1.000000
>>
>>
>>
>>
>>
>> BoomerAMG SOLVER PARAMETERS:
>>
>>   Maximum number of cycles:         6
>>   Stopping Tolerance:               1.000000e-05
>>   Cycle type (1 = V, 2 = W, etc.):  1
>>
>>   Relaxation Parameters:
>>    Visiting Grid:                     down   up  coarse
>>    Visiting Grid:                     down   up  coarse
>>             Number of partial sweeps:   20   20     1
>>    Type 0=Jac, 3=hGS, 6=hSGS, 9=GE:      6    6     6
>>    Point types, partial sweeps (1=C, -1=F):
>>                   Pre-CG relaxation (down):   1  -1   1  -1   1
>> -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1
>>                    Post-CG relaxation (up):  -1   1  -1   1  -1
>> 1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1
>>                              Coarsest grid:   0
>>
>>  Relaxation Weight 0.600000 level 0
>>  Outer relaxation weight 0.600000 level 0
>>  Output flag (print_level): 3
>>
>>
>> AMG SOLUTION INFO:
>>                                             relative
>>                residual        factor       residual
>>                --------        ------       --------
>>     Initial    7.914114e-03                 1.000000e+00
>>     Cycle  1   8.320104e-03    1.051300     1.051300e+00
>>     Cycle  2   7.767803e-03    0.933618     9.815126e-01
>>     Cycle  3   6.570752e-03    0.845896     8.302575e-01
>>     Cycle  4   5.836166e-03    0.888204     7.374377e-01
>>     Cycle  5   6.593214e-03    1.129717     8.330956e-01
>>     Cycle  6   8.191117e-03    1.242356     1.035001e+00
>>
>>
>> ==============================================
>>  NOTE: Convergence tolerance was not achieved
>>       within the allowed 6 V-cycles
>> ==============================================
>>
>>  Average Convergence Factor = 1.005750
>>
>>      Complexity:    grid = 1.000000
>>                 operator = 1.000000
>>                    cycle = 1.000000
>>
>>
>>
>>
>>
>> BoomerAMG SOLVER PARAMETERS:
>>
>>   Maximum number of cycles:         6
>>   Stopping Tolerance:               1.000000e-05
>>   Cycle type (1 = V, 2 = W, etc.):  1
>>
>>   Relaxation Parameters:
>>    Visiting Grid:                     down   up  coarse
>>    Visiting Grid:                     down   up  coarse
>>             Number of partial sweeps:   20   20     1
>>    Type 0=Jac, 3=hGS, 6=hSGS, 9=GE:      6    6     6
>>    Point types, partial sweeps (1=C, -1=F):
>>                   Pre-CG relaxation (down):   1  -1   1  -1   1
>> -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1
>>                    Post-CG relaxation (up):  -1   1  -1   1  -1
>> 1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1
>>                              Coarsest grid:   0
>>
>>  Relaxation Weight 0.600000 level 0
>>  Outer relaxation weight 0.600000 level 0
>>  Output flag (print_level): 3
>>
>>
>> AMG SOLUTION INFO:
>>                                             relative
>>                residual        factor       residual
>>                --------        ------       --------
>>     Initial    5.553016e-03                 1.000000e+00
>>     Cycle  1   6.417482e-03    1.155675     1.155675e+00
>>     Cycle  2   7.101926e-03    1.106653     1.278931e+00
>>     Cycle  3   7.077471e-03    0.996556     1.274527e+00
>>     Cycle  4   6.382835e-03    0.901853     1.149436e+00
>>     Cycle  5   5.392392e-03    0.844827     9.710744e-01
>>     Cycle  6   4.674173e-03    0.866809     8.417359e-01
>>
>>
>> ==============================================
>>  NOTE: Convergence tolerance was not achieved
>>       within the allowed 6 V-cycles
>> ==============================================
>>
>>  Average Convergence Factor = 0.971694
>>
>>      Complexity:    grid = 1.000000
>>                 operator = 1.000000
>>                    cycle = 1.000000
>>
>>
>>
>>
>>
>> BoomerAMG SOLVER PARAMETERS:
>>
>>   Maximum number of cycles:         6
>>   Stopping Tolerance:               1.000000e-05
>>   Cycle type (1 = V, 2 = W, etc.):  1
>>
>>   Relaxation Parameters:
>>    Visiting Grid:                     down   up  coarse
>>    Visiting Grid:                     down   up  coarse
>>             Number of partial sweeps:   20   20     1
>>    Type 0=Jac, 3=hGS, 6=hSGS, 9=GE:      6    6     6
>>    Point types, partial sweeps (1=C, -1=F):
>>                   Pre-CG relaxation (down):   1  -1   1  -1   1
>> -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1
>> -1   1  -1   1  -1
>>                    Post-CG relaxation (up):  -1   1  -1   1  -1
>> 1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1  -1   1  -1   1  -1   1  -1   1  -1
>> 1  -1   1  -1   1
>>                              Coarsest grid:   0
>>
>>  Relaxation Weight 0.600000 level 0
>>  Outer relaxation weight 0.600000 level 0
>>  Output flag (print_level): 3
>>
>>
>> AMG SOLUTION INFO:
>>                                             relative
>>                residual        factor       residual
>>                --------        ------       --------
>>     Initial    3.846663e-03                 1.000000e+00
>>     Cycle  1   4.136662e-03    1.075390     1.075390e+00
>>     Cycle  2   4.463861e-03    1.079097     1.160450e+00
>>     Cycle  3   4.302262e-03    0.963798     1.118440e+00
>>     Cycle  4   4.175328e-03    0.970496     1.085441e+00
>>     Cycle  5   4.735871e-03    1.134251     1.231163e+00
>>     Cycle  6   5.809054e-03    1.226607     1.510154e+00
>>
>>
>> ==============================================
>>  NOTE: Convergence tolerance was not achieved
>>       within the allowed 6 V-cycles
>> ==============================================
>>
>>  Average Convergence Factor = 1.071117
>>
>>      Complexity:    grid = 1.000000
>>                 operator = 1.000000
>>                    cycle = 1.000000
>>
>>
>>>
>>>   Nicolas,
>>>
>>>    I have added support for both 1 and/or 2 to petsc-dev
>>> http://www-unix.mcs.anl.gov/petsc/petsc-as/developers/index.html
>>>
>>> The two new options are -pc_hypre_boomeramg_nodal_coarsen and
>>> -pc_hypre_boomeramg_nodal_relaxation <optional:n) where the optional
>>> argument n indicates the levels which the SmoothNumLevels() sets.
>>>
>>> I have not tested the code so please let me know what problems you
>>> have.
>>>
>>>   Allison,
>>>
>>>     Thank you very much for the clarifications.
>>>
>>>    Barry
>>>
>>> On Wed, 25 Apr 2007, Allison Baker wrote:
>>>
>>>> Hi Barry and Nicolas,
>>>>
>>>> To clarify,
>>>>
>>>> HYPRE_BoomerAMGSetNumFunctions(solver, int num_functions) tells
>>>> AMG to
>>>> solve a
>>>> system of equations with the specified number of functions/
>>>> unknowns. The
>>>> default AMG scheme to solve a PDE system is the "unknown"
>>>> approach. (The
>>>> coarsening and interpolation are determined by looking at each
>>>> unknown/function independently - therefore you can imagine that the
>>>> coarse
>>>> grids are generally not the same for each variable. This approach
>>>> generally
>>>> works well unless you have strong coupling between unknowns.)
>>>>
>>>> HYPRE_BOOMERAMGSetNodal(solver, int nodal ) tells AMG to coarsen
>>>> such
>>>> that
>>>> each variable has the same coarse grid - sometimes this is more
>>>> "physical" for
>>>> a particular problem. The value chosen here for nodal determines how
>>>> strength
>>>> of connection is determined between the coupled system.  I suggest
>>>> setting
>>>> nodal = 1, which uses a Frobenius norm.  This does NOT tell AMG
>>>> to use
>>>> nodal
>>>> relaxation.
>>>>
>>>> If you want to use nodal relaxation in hypre there are two choices:
>>>>
>>>> (1) If you call  HYPRE_BOOMERAMGSetNodal, then you can
>>>> additionally do
>>>> nodal
>>>> relaxation via the schwarz smoother option in hypre.  I did not
>>>> implement this
>>>> in the Petsc interface, but it could be done easy enough. The
>>>> following
>>>> four
>>>> functions need to be called:
>>>>
>>>>      HYPRE_BoomerAMGSetSmoothType(solver, 6);
>>>>      HYPRE_BoomerAMGSetDomainType(solver, 1);
>>>>      HYPRE_BoomerAMGSetOverlap(solver, 0);
>>>>      HYPRE_BoomerAMGSetSmoothNumLevels(solver, num_levels);  (Set
>>>> num_levels
>>>> to number of levels on which you  want nodal smoothing, i.e.
>>>> 1=just the
>>>> fine
>>>> grid, 2= fine grid and the grid below, etc.  I find that doing nodal
>>>> relaxation on just the finest level is generally sufficient.)
>>>> Note that
>>>> the
>>>> interpolation scheme used will be the same as in the unknown
>>>> approach -
>>>> so
>>>> this is what we call a hybrid systems method.
>>>>
>>>> (2) You can do both nodal smoothing and a nodal interpolation
>>>> scheme.
>>>> While
>>>> this is currently implemented in 2.0.0, it is not advertised (i.e.,
>>>> mentioned
>>>> at all in the user's manual) because it is not yet implemented very
>>>> efficiently (the fine grid matrix is converted to a block matrix
>>>> - and
>>>> both
>>>> are stored), and we have not found it to be as effective as
>>>> advertised
>>>> elsewhere (this is an area of current research for us)..... If
>>>> you want
>>>> to try
>>>> it anyway, let me know and I will provide more info.
>>>>
>>>> Hope this helps,
>>>> Allison
>>>>
>>>>
>>>> Barry Smith wrote:
>>>>>   Nicolas,
>>>>>
>>>>> On Wed, 25 Apr 2007, Nicolas Bathfield wrote:
>>>>>
>>>>>
>>>>>> Dear Barry,
>>>>>>
>>>>>> Using MatSetBlockSize(A,5) improved my results greatly. Boomemramg
>>>> is now
>>>>>> solving the system of equations.
>>>>>>
>>>>>
>>>>>    Good
>>>>>
>>>>>
>>>>>> Still, the equations I solve are coupled, and my discretization
>>>> scheme is
>>>>>> meant for a non-segregated solver. As a consequence (I believe),
>>>> boomeramg
>>>>>> still diverges.
>>>>>>
>>>>>
>>>>>   How can "Boomeramg be now solving the system of equations" but
>>>>> also
>>>>> diverge? I am so confused.
>>>>>
>>>>>
>>>>>>  I would therefore like to use the nodal relaxation in
>>>>>> boomeramg (the hypre command is HYPRE_BOOMERAMGSetNodal) in
>>>>>> order to
>>>>>> couple the coarse grid choice for all my variables.
>>>>>>
>>>>>
>>>>>   I can add this this afternoon.
>>>>>
>>>>>   I have to admit I do not understand the difference between
>>>>> HYPRE_BOOMERAMGSetNodal() and hypre_BoomerAMGSetNumFunctions(). Do
>>>> you?
>>>>>
>>>>>    Barry
>>>>>
>>>>>> How can I achieve this from PETSc?
>>>>>>
>>>>>> Best regards,
>>>>>>
>>>>>> Nicolas
>>>>>>
>>>>>>
>>>>>>>   From PETSc MPIAIJ matrices you need to set the block size of
>>>>>>> the
>>>>>>> matrix
>>>>>>> with MatSetBlockSize(A,5) after you have called MatSetType() or
>>>>>>> MatCreateMPIAIJ(). Then HYPRE_BoomerAMGSetNumFunctions() is
>>>>>>> automatically
>>>>>>> called by PETSc.
>>>>>>>
>>>>>>>    Barry
>>>>>>>
>>>>>>> The reason this is done this way instead of as
>>>>>>> -pc_hypre_boomeramg_block_size <bs> is the idea that hypre will
>>>> use the
>>>>>>> properties of the matrix it is given in building the
>>>> preconditioner so
>>>>>>> the user does not have to pass those properties in seperately
>>>> directly
>>>>>>> to hypre.
>>>>>>>
>>>>>>>
>>>>>>> On Fri, 13 Apr 2007, Shaman Mahmoudi wrote:
>>>>>>>
>>>>>>>
>>>>>>>> Hi,
>>>>>>>>
>>>>>>>> int HYPRE_BoomerAMGSetNumFunctions (.....)
>>>>>>>>
>>>>>>>> sets the size of the system of PDEs.
>>>>>>>>
>>>>>>>> With best regards, Shaman Mahmoudi
>>>>>>>>
>>>>>>>> On Apr 13, 2007, at 2:04 PM, Shaman Mahmoudi wrote:
>>>>>>>>
>>>>>>>>
>>>>>>>>> Hi Nicolas,
>>>>>>>>>
>>>>>>>>> You are right. hypre has changed a lot since the version I
>>>> used.
>>>>>>>>>
>>>>>>>>> I found this interesting information:
>>>>>>>>>
>>>>>>>>> int HYPRE_BOOMERAMGSetNodal(....)
>>>>>>>>>
>>>>>>>>> Sets whether to use the nodal systems version. Default is 0.
>>>>>>>>>
>>>>>>>>> Then information about smoothers:
>>>>>>>>>
>>>>>>>>> One interesting thing there is this,
>>>>>>>>>
>>>>>>>>> HYPRE_BoomerAMGSetDomainType(....)
>>>>>>>>>
>>>>>>>>> 0 - each point is a domain (default)
>>>>>>>>> 1 each node is a domain (only of interest in systems AMG)
>>>>>>>>> 2 ....
>>>>>>>>>
>>>>>>>>> I could not find how you define the nodal displacement
>>>> ordering. But
>>>>>>>>>
>>>>>>>> it
>>>>>>>>
>>>>>>>>> should be there somewhere.
>>>>>>>>>
>>>>>>>>> I read the reference manual for hypre 2.0
>>>>>>>>>
>>>>>>>>> With best regards, Shaman Mahmoudi
>>>>>>>>>
>>>>>>>>>
>>>>>>>>> On Apr 13, 2007, at 1:40 PM, Nicolas Bathfield wrote:
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>> Dear Shaman,
>>>>>>>>>>
>>>>>>>>>> As far as I could understand, there is a BoomerAMG?s systems
>>>> AMG
>>>>>>>>>>
>>>>>>>> version
>>>>>>>>
>>>>>>>>>> available. This seems to be exactly what I am looking for,
>>>> but I
>>>>>>>>>>
>>>>>>>> just
>>>>>>>>
>>>>>>>>>> don't know how to access it, either through PETSc or
>>>> directly.
>>>>>>>>>>
>>>>>>>>>> Best regards,
>>>>>>>>>>
>>>>>>>>>> Nicolas
>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>>> Hi,
>>>>>>>>>>>
>>>>>>>>>>> You want to exploit the structure of the model?
>>>>>>>>>>> As far as I know, boomeramg can not treat a set of rows or
>>>>>>>>>>> blocks
>>>>>>>>>>>
>>>>>>>> as
>>>>>>>>
>>>>>>>>>>> a molecule,  a so called block-smoother?
>>>>>>>>>>> ML 2.0 smoothed aggregation does support it.
>>>>>>>>>>>
>>>>>>>>>>> With best regards, Shaman Mahmoudi
>>>>>>>>>>>
>>>>>>>>>>> On Apr 13, 2007, at 10:45 AM, Nicolas Bathfield wrote:
>>>>>>>>>>>
>>>>>>>>>>>
>>>>>>>>>>>> Hi,
>>>>>>>>>>>>
>>>>>>>>>>>> I am solving the Navier-stokes equations and try to use
>>>> Hypre
>>>>>>>>>>>> as
>>>>>>>>>>>> preconditioner.
>>>>>>>>>>>> Until now, I used PETSc as non-segregated solver and it
>>>> worked
>>>>>>>>>>>> perfectly.
>>>>>>>>>>>> Things got worse when I decided to use Boomeramg
>>>> (Hypre).
>>>>>>>>>>>> As I solve a system of PDEs,  each cell is represented
>>>> by 5
>>>>>>>>>>>> rows
>>>>>>>>>>>>
>>>>>>>> in my
>>>>>>>>
>>>>>>>>>>>> matrix (I solve for 5 variables). PETSc handles that
>>>> without
>>>>>>>>>>>>
>>>>>>>> problem
>>>>>>>>
>>>>>>>>>>>> apparently, but the coarsening scheme of Boomeramg needs
>>>> more
>>>>>>>>>>>>
>>>>>>>> input in
>>>>>>>>
>>>>>>>>>>>> order to work properly.
>>>>>>>>>>>>
>>>>>>>>>>>> Is there an option in PETSc to tell HYPRE that we are
>>>> dealing
>>>>>>>>>>>>
>>>>>>>> with a
>>>>>>>>
>>>>>>>>>>>> system of PDEs? (something like:
>>>> -pc_hypre_boomeramg_...)
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>> Thanks for your help.
>>>>>>>>>>>>
>>>>>>>>>>>> Best regards,
>>>>>>>>>>>>
>>>>>>>>>>>> Nicolas
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>> --
>>>>>>>>>>>> Nicolas BATHFIELD
>>>>>>>>>>>> Chalmers University of Technology
>>>>>>>>>>>> Shipping and Marine Technology
>>>>>>>>>>>> phone: +46 (0)31 772 1476
>>>>>>>>>>>> fax:   +46 (0)31 772 3699
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>
>>>>>>>>>> --
>>>>>>>>>> Nicolas BATHFIELD
>>>>>>>>>> Chalmers University of Technology
>>>>>>>>>> Shipping and Marine Technology
>>>>>>>>>> phone: +46 (0)31 772 1476
>>>>>>>>>> fax:   +46 (0)31 772 3699
>>>>>>>>>>
>>>>>>>>>>
>>>>>>
>>>>>>
>>>>>
>>>>>
>>>>
>>>>
>>>
>>>
>>
>>
>> --
>> Nicolas BATHFIELD
>> Chalmers University of Technology
>> Shipping and Marine Technology
>> phone: +46 (0)31 772 1476
>> fax:   +46 (0)31 772 3699
>>
>
>


-- 
Nicolas BATHFIELD
Chalmers University of Technology
Shipping and Marine Technology
phone: +46 (0)31 772 1476
fax:   +46 (0)31 772 3699




More information about the petsc-users mailing list