[Trilinos-Users] [EXTERNAL] Re: Xpetra wrappers around existing Tpetra objects
Martin Vymazal
martin.vymazal at vki.ac.be
Fri Jun 3 14:42:24 EDT 2016
Hi Andrey,
yes, that's exactly what I meant. Thank you.
Best regards,
Martin
On 2016-06-03 18:22, Andrey Prokopenko wrote:
> Hi Martin,
>
> If "MueLu as a solver' you mean a standalone solver, then you just
> need to add the following lines to what Jonathan wrote:
>
> Tpetra::MultiVector<> rhs, initial_guess;
> int max_its = 40;
>
> M->IsPreconditioner(false);
> M->Iterate(rhs, initial_guess, max_its);
>
> -Andrey
>
> On 06/03/2016 10:43 AM, Martin Vymazal wrote:
>> Hello Jonathan,
>>
>> thank you very much, this helps. I think that's all I need, but I was
>> wondering about one more thing - is there a similarly easy way to use
>> MueLu as a solver?
>>
>> Best regards,
>>
>> Martin
>>
>>
>> On 2016-06-03 17:28, Jonathan Hu wrote:
>>> Hi Martin,
>>>
>>> Yes, there are easier ways to use MueLu with an existing
>>> Tpetra::CrsMatrix. You do not need to convert to Xpetra. Something
>>> like this will work in your case:
>>>
>>> Teuchos::RCP<Tpetra::CrsMatrix<> > A;
>>> ParameterList paramList;
>>> //set the MueLu options
>>> Teuchos::RCP<MueLu::TpetraOperator> M =
>>> MueLu::CreateTpetraPreconditioner(A, paramList);
>>>
>>> M is a Tpetra::Operator and can be used directly with Belos.
>>>
>>> This is described in more detail in the user's guide
>>> (muelu/doc/UsersGuide/mueluguide.pdf) on page 22. Thank you for the
>>> feedback on the documentation ... I'll see if we can make the
>>> migration path clearer for first-time users.
>>>
>>> Regards,
>>> Jonathan
>>>
>>> Martin Vymazal wrote on 06/03/2016 08:10 AM:
>>>> Hello,
>>>>
>>>> I have an RCP pointer to Tpetra::CrsMatrix called A and I would
>>>> like to pass this to MueLu Hierarchy object. It seems that the only
>>>> option is to wrap A, the right-hand side vector and the vector of
>>>> unknows in Xpetra RCP pointers. The following compiles (I found a
>>>> similar code in the MueLU tutorial):
>>>>
>>>> Teuchos::RCP<Xpetra::CrsMatrix<>> exA = Teuchos::rcp(new
>>>> Xpetra::TpetraCrsMatrix<>(A));
>>>> Teuchos::RCP<Xpetra::CrsMatrixWrap<>> exAWrap = Teuchos::rcp(new
>>>> Xpetra::CrsMatrixWrap<>(exA));
>>>> RCP<Xpetra::Matrix<>> xpetra_A =
>>>> Teuchos::rcp_dynamic_cast<Xpetra::Matrix<>>(exAWrap);
>>>>
>>>> MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal> H(xpetra_A);
>>>> H.SetDefaultVerbLevel(MueLu::Medium);
>>>>
>>>> Could you please confirm that the first line only creates a shallow
>>>> copy of A? Is there any easier/more correct/more elegant way of
>>>> constructing xpetra_A? I'm also stuck at this point, because now I
>>>> need Xpetra-style pointers to the map that was used to construct
>>>> (Tpetra) A in order to construct Xpetra-style right-hand side
>>>> multivector (or pointer to it). How do I get Xpetra vectors from
>>>> Tpetra vectors, given the fact that the map of A is a Tpetra map?
>>>>
>>>> Thank you for your help.
>>>>
>>>> Martin
>>>>
>>>> P. S. MueLu looks interesting, but the documentation (especially the
>>>> Xpetra part) could be a bit more extensive ... It looks like there's
>>>> a lot of attention paid to users who migrate from Epetra to Tpetra,
>>>> but what about those, who already have codebase which is exclusively
>>>> Tpetra-based? Are they supposed to move everything to Xpetra in
>>>> order to be able to use MueLu?
>>>> _______________________________________________
>>>> Trilinos-Users mailing list
>>>> Trilinos-Users at trilinos.org
>>>> https://trilinos.org/mailman/listinfo/trilinos-users
>>>>
>>>
>>> _______________________________________________
>>> Trilinos-Users mailing list
>>> Trilinos-Users at trilinos.org
>>> https://trilinos.org/mailman/listinfo/trilinos-users
>> _______________________________________________
>> Trilinos-Users mailing list
>> Trilinos-Users at trilinos.org
>> https://trilinos.org/mailman/listinfo/trilinos-users
More information about the Trilinos-Users
mailing list