Hi again Feilpe,

Regarding the lack of documentation, I’ve opened an issue on Github to track 
this. You can find that here: https://github.com/dealii/dealii/issues/12466 
<https://github.com/dealii/dealii/issues/12466>

Best,
Jean-Paul

> On 15. Jun 2021, at 12:57, Jean-Paul Pelteret <jppelte...@gmail.com> wrote:
> 
> Hi Feilpe,
> 
> Firstly, I agree that the documentation is very light on details on how to 
> use the linear operators with Trilinos linear algebra types. We can 
> definitely improve this, and any contributions to that effect would be 
> greatly appreciated!
> 
> Right now, I can direct you to a few of the tests that use Trilinos LA in 
> conjunction with the inverse operator, so that you can compare what you have 
> and figure out what the problematic differences are. There is this one, for 
> instance
> https://github.com/dealii/dealii/blob/master/tests/lac/linear_operator_12a.cc#L323-L341
>  
> <https://github.com/dealii/dealii/blob/master/tests/lac/linear_operator_12a.cc#L323-L341>
> that looks like a similar setup to yours, and
> https://github.com/dealii/dealii/blob/master/tests/lac/schur_complement_04.cc#L126-L137
>  
> <https://github.com/dealii/dealii/blob/master/tests/lac/schur_complement_04.cc#L126-L137>
> that uses a Trilinos::SolverCG (as opposed to deal.II’s solver).
> 
> Comparing to both of these, I think that the important point might be that 
> the preconditioner must be declared (out-of-line) before the 
> inverse_operation() function is called. The linear operators typically expect 
> the lifetime of the LA objects to exceed that of the operators, and so you 
> have to first create the matrix or preconditioner and then pass it to the 
> linear operators. This is similar to what you’d done when setting up op_M, 
> for example. The case of passing in a deal::PreconditionerIdentity() for 
> serial operations is a special case, and I don’t think that we’ve duplicated 
> that for TrilinosWrappers:: PreconditionerIdentity(). Maybe that could be 
> improved too.
> 
> I hope that with this single change you’d be able to get your program 
> running. If not, then please do let us know so that we can try to help 
> further.
> 
> Best,
> Jean-Paul
> 
> 
>> On 14. Jun 2021, at 19:18, Juan Felipe Giraldo <jfgiral...@gmail.com 
>> <mailto:jfgiral...@gmail.com>> wrote:
>> 
>> Hello everyone,
>> 
>> I have implemented a residual-minimization framework that somehow is similar 
>> to DPG. I want to extend my results by using parallelization using MPI with 
>> PETSc or Trilinos. 
>> So far, I have solved the saddle point problem using the Schur complement 
>> exactly how it is described in step 20. Now, I am trying to replicate 
>> exactly the same solver but using the MPI wrappers and the linear operators.
>> 
>> The problem is that when I am trying to implement the inverse_operator to 
>> compute the Preconditioner of the Schur complement, I get an error saying 
>> that the functions are not matching "inverse_operator(op_aS, solver_aS, 
>> TrilinosWrappers::PreconditionIdentity())."
>> 
>> There is no much documentation about linear operators in parallel solvers, 
>> so if anyone has any suggestion on how to fix this problem, it would be well 
>> appreciated. 
>> 
>> I have pasted the complete function in below:
>> 
>> 
>> template <int dim>
>> void FEMwDG<dim>::
>>   solve()
>>   {
>>     TimerOutput::Scope t(computing_timer, "solve");
>> 
>>     LA::MPI::PreconditionJacobi prec_M;
>>      LA::MPI::PreconditionJacobi::AdditionalData data;
>>      prec_M.initialize(system_matrix.block(0, 0), data);
>>     }
>> 
>>     auto &E = solution.block(0);
>>     auto &U = solution.block(1);
>>     const auto &L = system_rhs.block(0);
>>     const auto  M = linear_operator<LA::MPI::Vector>(system_matrix.block(0, 
>> 0));
>> 
>>     const auto op_M  = linear_operator(M, prec_M);
>>     const auto op_B   = 
>> linear_operator<LA::MPI::Vector>(system_matrix.block(0, 1));
>>     const auto op_BT = 
>> linear_operator<LA::MPI::Vector>(system_matrix.block(1, 0));
>> 
>>     ReductionControl          inner_solver_control2(5000,
>>                                                     1e-18 * 
>> system_rhs.l2_norm(),
>>                                                     1.e-2);
>> 
>>     SolverCG<LA::MPI::Vector> cg(inner_solver_control2);
>>     const auto op_M_inv = inverse_operator(M, cg, prec_M);
>> 
>>     const auto op_S = op_BT * op_M_inv * op_B;
>>     const auto op_aS = op_BT * linear_operator<LA::MPI::Vector>(prec_M) * 
>> op_B;
>> 
>>     IterationNumberControl   iteration_number_control_aS(30, 1.e-18);
>>     SolverCG<LA::MPI::Vector> solver_aS(iteration_number_control_aS);
>> 
>>     const auto preconditioner_S =
>>     inverse_operator(op_aS, solver_aS, 
>> TrilinosWrappers::PreconditionIdentity());
>>     const auto schur_rhs = op_BT * op_M_inv * L ;
>> 
>>     SolverControl            solver_control_S(2000, 1.e-12);
>>     SolverCG<LA::MPI::Vector> solver_S(solver_control_S);
>> 
>>     const auto op_S_inv = inverse_operator(op_S, solver_S, preconditioner_S 
>> );
>> 
>>     U = op_S_inv * schur_rhs;
>>     std::cout << solver_control_S.last_step()
>>               << " CG Schur complement iterations to obtain convergence."
>>               << std::endl;
>>     E = op_M_inv * (L - op_B * U);
>>      }
>> 
>> Thank you so much, 
>> 
>> Regards,
>> Felipe
>> 
>> 
>> 
>> 
>> 
>> -- 
>> The deal.II project is located at http://www.dealii.org/ 
>> <http://www.dealii.org/>
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en 
>> <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 
>> <mailto:dealii+unsubscr...@googlegroups.com>.
>> To view this discussion on the web visit 
>> https://groups.google.com/d/msgid/dealii/40b8db07-5243-4603-825b-9e7850580206n%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/40b8db07-5243-4603-825b-9e7850580206n%40googlegroups.com?utm_medium=email&utm_source=footer>.
> 

-- 
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.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/B0189A0A-376F-4A9F-A1ED-DA6743B031A8%40gmail.com.

Reply via email to