[Trilinos-Users] ML package (standalone) not scaling

Deb Rajdeep debr at student.ethz.ch
Tue Jul 19 11:11:44 EDT 2016


Hello Sir/Madam,
I am a user of ML package of Trilinos. I am using it for an iterative matrix solver problem. The matrix has 9 point stencil but diagonally dominant structure. The serial code is working fine and converging. But I have issues with parallelization.

My problem is that even though parallel code is converging but its never scaling beyond 4 or 5 processors. So following is the description of my job.

The domain size is 320 by 320 grid nodes. If I use four processors, then each processor is having a matrix row size of 40 by 40 grid nodes. Which means 1600 row sizes. Then there are ghost cells whose information is transferred. One of the snapshot of converging solution with 4 processors is


ML_Smoother_NewGS: can't get Crs data pointers (return code -1, mat_type = 0)
Iter       ||res_i||_2**    ||res_i||/||res_i+1||
   1        1.501e-03
   3        5.931e-05            2.898e-01
   5        3.834e-06            2.510e-01
   7        2.588e-07            2.662e-01
   9        2.722e-08            3.528e-01

**Residual norm taken after multigrid pre_smoothing step.


It takes much less than time than 2 processors to compute. However if I use 8 or any higher number of cores than it takes more time to compute. Can you please guide me about if I can resolve this issue by changing of some parameters or what can possibly be the issue.

Following are the parameters set up before running ML_Iterate()


   N_grids = 20;
   ML_Create         (&ml_object, N_grids);
   int proc = ml_object->comm->ML_mypid;

   ML_Init_Amatrix      (ml_object, 0, localRowSize, localRowSize,tmpSys);
   ML_Set_Amatrix_Getrow(ml_object, 0, Poisson_getrow, Poisson_comm, localColSize);
   ML_Set_Amatrix_Matvec(ml_object, 0, Poisson_matvec);

   ML_Set_PrintLevel(1);
   ML_Set_ResidualOutputFrequency(ml_object,iterationStep);
   ML_Set_MaxIterations(ml_object, iterations);
   ML_Set_Tolerance(ml_object, residual);

   ML_Aggregate_Create(&agg_object);
   ML_Aggregate_Set_MaxCoarseSize(agg_object,1);
   N_levels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object);

   ML_Gen_Smoother_GaussSeidel(ml_object, ML_ALL_LEVELS, ML_BOTH, 8, ML_DEFAULT);
   ML_Gen_Solver    (ml_object, ML_MGV, 0, N_levels);

   ML_Iterate(ml_object,sol,rhs);
   ML_Aggregate_Destroy(&agg_object);
   ML_Destroy(&ml_object);


Thanks and regards,
Rajdeep




-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20160719/24fdc141/attachment.html>


More information about the Trilinos-Users mailing list