Dear Jean-Paul,
Thank you, everything works. Right now for simple test cases I have managed 
to get interesting results:

   - Going from ILU to AMG leads to speed-up of 4X on the same number of 
   processor.
   -  More importantly, my number of GMRES iteration remains extremely low 
   (below 15) and is practically independent of the number of cells (10 
   iterations with 128x128, 12 with 256x256, etc.). If my multigrid solver 
   classes are not too far behind me, this should be the expected behavior

I need to test it more on more complex test cases, but once this is done, I 
will open a pull request and try to make an interesting tutorial out of it 
this summer.
I just wanted to say thanks for your amazing contribution, these small 
modifications really open up a lot of possibilities in terms of exploring 
the AMG.

With this, I have weeks of fun ahead :)! 

Here is a small final example for the setting of the ILU smoother and 
coarse parameters.

    TrilinosWrappers::PreconditionAMG preconditioner;

  std::vector<std::vector<bool> > constant_modes;
  std::vector<bool>  velocity_components (dim+1,true);
  DoFTools::extract_constant_modes (dof_handler, velocity_components,
                                    constant_modes);
  TrilinosWrappers::PreconditionAMG::AdditionalData amg_data;
  amg_data.constant_modes = constant_modes;

  const bool            elliptic=false;
  const bool            higher_order_elements = false;
  const unsigned int      n_cycles = 2;
  const bool            w_cycle = true;
  const double      aggregation_threshold = 1e-10;
  const unsigned int      smoother_sweeps = 2;
  const unsigned int      smoother_overlap = 1;
  const bool            output_details = false;
  const char *      smoother_type = "ILU";
  const char *      coarse_type = "ILU";
  TrilinosWrappers::PreconditionAMG::AdditionalData preconditionerOptions(
        elliptic,
        higher_order_elements,
        n_cycles,
        w_cycle,
        aggregation_threshold,
        constant_modes,
        smoother_sweeps,
        smoother_overlap,
        output_details,
        smoother_type,
        coarse_type
        );

  Teuchos::ParameterList            parameter_ml;
  std::unique_ptr< Epetra_MultiVector > distributed_constant_modes;
  preconditionerOptions.set_parameters(parameter_ml, 
distributed_constant_modes, system_matrix);
  const double ilu_fill=linearSolverParameters.ilut_fill;
  const double ilu_atol=linearSolverParameters.ilut_atol ;
  const double ilu_rtol=linearSolverParameters.ilut_rtol;
  parameter_ml.set("smoother: ifpack level-of-fill",ilu_fill);
  parameter_ml.set("smoother: ifpack absolute threshold",ilu_atol);
  parameter_ml.set("coarse: ifpack level-of-fill",ilu_fill);
  parameter_ml.set("coarse: ifpack absolute threshold",ilu_atol);
  preconditioner.initialize(system_matrix,parameter_ml);



On Tuesday, 2 April 2019 05:55:26 UTC-4, Jean-Paul Pelteret wrote:
>
> Dear Bruno,
>
> These are the basic steps that you need to take:
>
> 1. Create and initialise a PreconditionAMG::AdditionalData object 
> (“additional_data”) as per usual. Here you set up the constant modes as is 
> needed using DoFTools::extract_constant_modes().
> 2. Create a Teuchos::ParameterList (“parameter_ml”) and std::unique_ptr< 
> Epetra_MultiVector > (“distributed_constant_modes”). Neither of these need 
> be initialised with any data.
> 3. Now call additional_data.set_parameters(parameter_ml, 
> distributed_constant_modes, matrix); .  This initialises “ parameter_ml” 
> and “”distributed_constant_modes” for you, and transcribes the settings 
> stored in “additional_data” into the equivalent “parameters_ml”. The 
>  “distributed_constant_modes” is internally referenced by “ parameters_ml”, 
> so thats why you need that vector lying around.
> 4. Lastly, you initialise the preconditioner using ” parameter_ml” for the 
> “matrix” that you sent into “ additional_data.set_parameters()”, i.e. 
> something like this: const PreconditionAMG prec (matrix, parameter_ml); .
>
> As for the second point, you should be able to find a link to the ML 
> documentation in the description to this function:
>
> https://dealii.org/developer/doxygen/deal.II/classTrilinosWrappers_1_1PreconditionAMG.html#a9703cc100be147a47358cf99fc53beb0
>
> You can also search for this technical report which should outline most if 
> not all of the possible settings:
>
> *TechReport (Gee2006a)*
> Gee, M. W.; Siefert, C. M.; Hu, J. J.; Tuminaro, R. S. & Sala, M. G.
> ML 5.0 Smoothed Aggregation User's Guide 
> *Sandia National Laboratories, * *Sandia National Laboratories, * *2006* 
>
> Best,
> Jean-Paul
>
> On 01 Apr 2019, at 17:53, Bruno Blais <blais...@gmail.com <javascript:>> 
> wrote:
>
> Dear Jean-Paul,
> I've compiled the new DEALII version with your changes, but I am unsure of 
> how to use the changes correctly.
> If I understand correctly, you can now set Teuchos parameters after you 
> have initialized the TrilinosWrappers::PreconditionAMG.
> This is achieved by creating a Teuchos::ParameterList object and using it 
> in the set_parameters public function of the 
> Preconditioner::AdditionalParameters.
> So for instance in my case it would be something like:
>
>
> // Still need to figure out the syntax for the additional_parameters
>  Teuchos::ParameterList            additional_parameter("smoother: params"
> );
>  
>
> And then I would set the parameters in the preconditioner additional 
> parameters
>
> I have two things i am unsure:
> 1) For the set the additional parameters I need a 
> std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes . I presume 
> this vector is created in the initialization from the std::vector of 
> constant boolean modes? Is there a way to grab it again and reuse it? I do 
> not understand why I need to specify the constant modes for every change of 
> parameters.
> 2) Is there any documentation on the various parameter list. What I would 
> like to change is the atol of the ILU smoother and ILU preconditions. I 
> think this will be in the smoother: params, but I am unsure how to proceed 
> from there. I have found the MULUE documentation, but it seems to be more 
> tailored around their XML interface?
>
>
> Thank you so much,
> Bruno
>
>
>
> On Tuesday, 12 March 2019 05:06:03 UTC-4, Jean-Paul Pelteret wrote:
>>
>> Dear Bruno,
>>
>> I believe that my proposed changes are working now, so you could always 
>> pull that branch from my repository and use it immediately. That might be 
>> the most convenient way for you to move forward.
>>
>> Thanks for being willing to share your findings wrt. AMG parameters. I 
>> think that would be really useful. I also think that encapsulating a 
>> parameter study for NS + AMG/GMG within a tutorial could be interesting and 
>> valuable. Maybe you could open up an issue on the GitHub repository and we 
>> could all discuss the specifics there.
>>
>> Best,
>> Jean-Paul
>>
>> On 11 Mar 2019, at 14:25, Bruno Blais <blais...@gmail.com> wrote:
>>
>>  Dear Jean-Paul and Wolfgang,
>> Is it better if I try to go with the Teuchos way or should I wait until 
>> Jean-Paul finishes the merge on his branch ? I think the latter might be a 
>> good option.
>> In all cases, I would be more than glad to share my experience in 
>> parameters tuning for the AMG solver with ILU type of smoother and coarse 
>> solve.
>>
>> On a slightly related note. Hopefully I will have finished the basis of 
>> my GLS stabilized navier-stokes solver by June or so. It already works 
>> pretty well in parallel with GMRES + ILU, but I want to test more what can 
>> be done in terms of multigrid and AMG.
>> In a way, it is very reminescent of a similar work done in MOOSE, but 
>> using Trilinos, P4est and deal.II. The MOOSE paper is here:
>> https://www.sciencedirect.com/science/article/pii/S0965997817310591
>>
>> I am more than willing to share everything about the solver. Do you 
>> believe this could be something that could make an interesting Tutorial? I 
>> could wrap-it around the study of the flow around a cylinder or a backward 
>> facing step to have some validation. Verification would be done using MMS.
>> Best
>> Bruno
>>
>> -- 
>> The deal.II project is located at http://www.dealii.org/
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en
>> --- 
>> You received this message because you are subscribed to the Google Groups 
>> "deal.II User Group" group.
>> To unsubscribe from this group and stop receiving emails from it, send an 
>> email to dea...@googlegroups.com.
>> For more options, visit https://groups.google.com/d/optout.
>>
>>
>>
> -- 
> The deal.II project is located at http://www.dealii.org/
> For mailing list/forum options, see 
> https://groups.google.com/d/forum/dealii?hl=en
> --- 
> You received this message because you are subscribed to the Google Groups 
> "deal.II User Group" group.
> To unsubscribe from this group and stop receiving emails from it, send an 
> email to dea...@googlegroups.com <javascript:>.
> For more options, visit https://groups.google.com/d/optout.
>
>
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to