[Trilinos-Users] Can we expect order-optimal behaviour from ML?

Ray Tuminaro rstumin at sandia.gov
Fri Dec 22 11:39:11 MST 2006


gunnaran at simula.no wrote:
> Dear ML experts
>
> My question is simply, can we expect order-optimal behavior for ML on
> elliptic and parabolic problems?
>
> My definition of order-optimality is: The cost of solving the problem
> grows linearly with the number of unknowns. This means that there has to
> be an upper bound of the number of iterations, regardless of the number of
> unknowns.
>
> Attached is two plots.
>
> cond_poisson_ml_P1_.eps is a simple poisson problem with linear finite
> elements on a unit square. The y-axis (kappa) is the condition number,
> estimated using conjugated gradient (described in the book by Saad) with
> ML as preconditioner. One V cycle is used, with symmetric Gauss-Seidel as
> point smoother. The condition number is oscillation, but I guess this is
> to be expected since not all grids can be equally well coarsened. The
> lower bound seems to be close to order-optimal, but the upper bound is far
> from it.
>
> cond_heat_ml_P1_s1_.eps is a simple heat problem on the unit square,
> discretized using linear finite elements and implicit Euler. The timestep
> is 0.1. Again, conjugated gradient is used to estimate the condition
> number. As we can see from the plot, we obviously have a O(log n) behavior
> with respect to the condition number, which gives us a O(n log n) behavior
> when it comes to computational cost.
>
> Geometric multigrid is known to be order-optimal, and it can be shown that
> for some grids, algebraic multigrid will have the same behavior. Should I
> expect order-optimal behavior, or is this O(n log n) what one should
> expect?
> The aggregation parameters are:
>
>  ML_Set_PrintLevel(print_level);
>  ML_Create(&ml, 20);
>  ML_Init_Amatrix       (ml, 0, n, n, (void *)csr_data);
>  ML_Set_Amatrix_Getrow (ml, 0, CSR_getrow, NULL, n);
>  ML_Set_Amatrix_Matvec (ml, 0, CSR_matvec);
>    cout << "Printing diagonal"  << endl;
>    for (int i = 0; i<n; i++) cout << A->diagonal[i] << endl;
>  if (not A->diagonal_computed)
>      A->computeDiag();
>  ML_Set_Amatrix_Diag   (ml, 0, n, A->diagonal);
>
>  ML_Aggregate_Create(&agg_object);
>  ML_Aggregate_Set_MaxCoarseSize(agg_object, 10);
>  N_levels = ML_Gen_MGHierarchy_UsingAggregation(ml, 0, ML_INCREASING,
> agg_object);
>  ML_Gen_Smoother_SymGaussSeidel(ml, ML_ALL_LEVELS, ML_BOTH, 1, 1.0);
>  ML_Gen_Solver    (ml, ML_MGV, 0, N_levels-1);
>
>
> The preconditioning is called using:
>
>     ML_Cycle_MG( &(ml->SingleLevel[ml->ML_finest_level]), d, c,
>             ML_NONZERO,ml->comm, ML_COMPUTE_RES_NORM,ml);
>
>
> Am I missing something important?
>
> Best regards
>
> Gunnar Staff
> ------------------------------------------------------------------------
>
> _______________________________________________
> Trilinos-Users mailing list
> Trilinos-Users at software.sandia.gov
> http://software.sandia.gov/mailman/listinfo/trilinos-users
Gunnar,

   The answer to your question is .... yes and no.

If I properly understand your plots, they are fairly typical of what we 
see.
That is, very modest growth in condition number as the mesh is refined. 
I haven't looked
at condition numbers recently, but we frequently see a slight increase 
in iterations
as the mesh is refined for Poisson problems.  In your case, it looks 
like the mesh is refined to something equivalent to (500 x 500) with a 
condition number
rising from about 1 to 4. Here is some data that I just collected for a 
finite difference Laplace operator:

mesh size         iters
=======    ========
 10x10               6
 30x30               7
 50x50               7
 70x70               7
 80x80               8
 90x90               8
100x100              7
120x120              8
140x140              8
160x160              7
180x180              8
280x280              7
380x380              9
480x480              8
580x580              8

There is a general trend toward iteration growth. We BELIEVE that this 
growth is
due to irregularities in the aggregates and could be rectified with 
better aggregation.
These aggregation irregularities are most likely the cause of the 
oscillations that you
see. Making perfect aggregates, however, is not easy. We are currently
working on variants which are a little less sensitive to aggregation ... 
but this
is on-going. There are a few routines that you can play with to try and get
different aggregates. There is an MIS aggregation, Metis aggregation, etc.
There is also something called:
    ML_Aggregate_Set_Phase3AggregateCreationAggressiveness
which will alter how aggregates are picked. Overall, however, you will
likely still see some modest iteration growth.

-Ray


-- 
Ray Tuminaro                      phone: (925) 294-2564
MS 9159                           fax:   (925) 294-2234
Sandia National Laboratories      email: rstumin at sandia.gov
PO Box 969                        http://www.cs.sandia.gov/~rstumin
Livermore, CA 94551



-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20061222/4af15baf/attachment.html


More information about the Trilinos-Users mailing list