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

gunnaran at simula.no gunnaran at simula.no
Tue Dec 19 00:50:44 MST 2006


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
-------------- next part --------------
A non-text attachment was scrubbed...
Name: cond_heat_ml_P1_s1_.eps
Type: application/postscript
Size: 148762 bytes
Desc: not available
Url : http://software.sandia.gov/mailman/private/trilinos-users/attachments/20061219/e683309a/cond_heat_ml_P1_s1_-0001.eps
-------------- next part --------------
A non-text attachment was scrubbed...
Name: cond_poisson_ml_P1_.eps
Type: application/postscript
Size: 149116 bytes
Desc: not available
Url : http://software.sandia.gov/mailman/private/trilinos-users/attachments/20061219/e683309a/cond_poisson_ml_P1_-0001.eps


More information about the Trilinos-Users mailing list