[Trilinos-Users] [EXTERNAL] ml with Anasazi
Thornquist, Heidi K
hkthorn at sandia.gov
Thu Jul 17 10:50:22 MDT 2014
Hi Sunghwan,
Anasazi's preconditioned eigensolvers, like LOBPCG and Block Davidson, can work directly with ML. For block Krylov-Schur, the ML preconditioner can be used with a spectral transformation, like shift-invert, to get the smallest eigenpairs.
Your example showing how the ML preconditioner is constructed is missing one crucial element that allows the Epetra_Operator generated by ML to work with Anasazi, this operator needs to be wrapped in an inverse operator so that calls to the Epetra_Operator method Apply actually call the ML preconditioner's ApplyInverse method. Without that, Anasazi will call the ML preconditioner's Apply method, which is not defined for this interface. See below.
Teuchos::RCP < Anasazi::BasicEigenproblem<double,MV,OP> > MyProblem =Teuchos::rcp( new Anasazi::BasicEigenproblem<double,MV,OP>(matrix, eigenvectors) );
Teuchos::RCP<ML_Epetra::MultiLevelPreconditioner> prec;
<
…
construct ML preconditioner, prec, as below
…
>
Teuchos::RCP<Epetra_Operator> PrecOp;
PrecOp = Teuchos::rcp( new Epetra_InvOperator(prec.get()) );
MyProblem->setPrec(PrecOp);
Thanks,
Heidi
--
Heidi K. Thornquist
Electrical Models & Simulation
Sandia National Laboratories
Albuquerque, NM 87185-1323
From: SungHwan Choi <sunghwanchoi91 at gmail.com<mailto:sunghwanchoi91 at gmail.com>>
Date: Thursday, July 17, 2014 7:15 AM
To: "trilinos-users at software.sandia.gov<mailto:trilinos-users at software.sandia.gov>" <trilinos-users at software.sandia.gov<mailto:trilinos-users at software.sandia.gov>>
Subject: [EXTERNAL] [Trilinos-Users] ml with Anasazi
Dear all,
Hi, I am an Anasazi user. I want to use a multilevel preconditioner. As far as, I know AztecOO work well with ML package. However, I am not sure multigrid like preconditioner can be applied on Anasazi.
The article titled "A comparison of eigensolvers for large-scale 3D modal analysis using AMG-preconditioned iterative methods" said that ML package can work with Anasazi solver but I couldn't find how ML can work with Anasazi.
Thus, applying my pre-knowledge, I coded it but It does not work well. Would you help me?
Teuchos::RCP < Anasazi::BasicEigenproblem<double,MV,OP> > MyProblem =Teuchos::rcp( new Anasazi::BasicEigenproblem<double,MV,OP>(matrix, eigenvectors) );
Teuchos::RCP<Epetra_Operator> PrecOp;
Teuchos::ParameterList MLList;
ML_Epetra::SetDefaults("DD",MLList);
MLList.set("smoother: pre or post", "post");
MLList.set("PDE equations", 1);
MLList.set("smoother: type","IFPACK");
MLList.set("smoother: ifpack type", "ILU");
MLList.set("smoother: ifpack overlap", 1);
MLList.sublist("smoother: ifpack list").set("fact: level-of-fill", 5);
ML_Epetra::MultiLevelPreconditioner* MLPrec =new ML_Epetra::MultiLevelPreconditioner(*matrix.get(), MLList, true);
PrecOp=Teuchos::rcp( MLPrec);
MyProblem->setPrec(PrecOp);
Sunghwan Choi
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://software.sandia.gov/pipermail/trilinos-users/attachments/20140717/5f4b749c/attachment.html>
More information about the Trilinos-Users
mailing list