[Trilinos-Users] Ifpack_Preconditioner

hkthorn hkthorn at sandia.gov
Thu Mar 16 15:22:42 MST 2006


Hi Bora,

I wasn't sure if you'd heard back from anyone about this problem, but here's
a response to what is going on:

Yes, you did what the Ifpack User's Guide says on pg 22-23.  However, upon
examination of the code, that's not going to compile.  That is because, as
your compilation error correctly states, the Inverse() method does not exist
for the abstract base class (Ifpack_Preconditioner).  It does exist for the
Ifpack_AdditiveSchwarz class.  However, when you declared the preconditioner

> Ifpack_Preconditioner *BJforSC;

you declared a pointer to an object of type Ifpack_Preconditioner.  To
access the Ifpack_AdditiveSchwarz specific functionality like Inverse()
you will have to dynamically cast the pointer, i.e.

> const Ifpack_ILUT* Inverse =
> dynamic_cast<Ifpack_AdditiveSchwarz<Ifpack_ILUT>* >(BJforSC)->Inverse();

You get back a const pointer, instead of a non-const pointer as the
documentation states.  Then you can continue to get the factors L() and U()
from Inverse.

This should solve your problem, but may I suggest:

1)  Declaring the preconditioner to be of type
Ifpack_AdditiveSchwarz<Ifpack_ILUT>.  This would remove the need for the
dynamic cast.  That is only good if you aren't trying to pass around an
all-purpose abstract object that will be generated later by some other part
of the code.  If that is the case, the dynamic cast that I suggested might
return a failure if the preconditioner is not of type
Ifpack_AdditiveSchwarz<Ifpack_ILUT>.

2)  Creating an object or a reference-counted pointer i.e.

Ifpack_AdditiveSchwarz<Ifpack_ILUT> BJforSC( SC, OverlapLevel );

Teuchos::RefCountPtr<Ifpack_AdditiveSchwarz<Ifpack_ILUT> > BJforSC =
 Teuchos::rcp( new Ifpack_AdditiveSchwarz<Ifpack_ILUT>(SC, OverlapLevel) );

This suggestion is merely for the sake of memory management :)

Hopefully this answers your question, if it hasn't been addressed already!

Thanks!

Heidi

-- 
Heidi K. Thornquist
Sandia National Laboratories
P.O. Box 5800 / MS 0316
Albuquerque, NM 87185
Office: (505)284-8426
Fax: (505)284-5451



On 3/14/06 8:41 PM, "Bora Ucar" <boraucar at gmail.com> wrote:

> Hi,
> 
> I am constructing ILUT preconditioner using the lines
> ------------------------------------------------------------------------------
> -----------------------------------------
> /*Matrix SC is already fillcomplete'd*/
> 
> OverlapLevel = 0;
> List.set("fact: drop tolerance", 1e-3);
> Ifpack_Preconditioner *BJforSC;
> BJforSC = new Ifpack_AdditiveSchwarz<Ifpack_ILUT>(SC, OverlapLevel);
> BJforSC->SetParameters(List);
> BJforSC->Initialize();
> BJforSC->Compute();
> 
> ------------------------------------------------------------------------------
> -----------------------------------------
> 
> I need to access to number of nonzeros in the approximate factors L and U.
> 
> Ifpack manual (SAND2005-0662) says the lines
> 
> Ifpack_ILUT* Inverse = BJforSC->Inverse();
> const Epetra_CrsMatrix& L = Inverse->L();
> 
> should work.
> 
> But, I got the following compilation error
> ------------------------------------------------------------------------------
> -----------------------------------------
> proPre.cpp: In constructor `proPre::proPre(Epetra_CrsMatrix&,
>    Epetra_CrsMatrix&, const Epetra_Comm&, const Epetra_Map&)':
> proPre.cpp:71: no matching function for call to
> `Ifpack_Preconditioner::Inverse
>    ()'
> ------------------------------------------------------------------------------
> -----------------------------------------
> 
> The file proPre.cpp: includes
> 
> #include "Ifpack_config.h"
> #include "Ifpack_ConfigDefs.h"
> #include "Ifpack.h"
> #include "Ifpack_AdditiveSchwarz.h"
> #include "Ifpack_ILUT.h"
> #include "Ifpack_Preconditioner.h"
> 
> Can you see any problems?
> 
> 
> _______________________________________________
> Trilinos-Users mailing list
> Trilinos-Users at software.sandia.gov
> http://software.sandia.gov/mailman/listinfo/trilinos-users



More information about the Trilinos-Users mailing list