[Trilinos-Users] ilut setting

Raymond Tuminaro rstumin at sandia.gov
Thu Sep 26 15:25:44 MDT 2013


I should probably point out that there is an unadvertised 'AZ_fixed_pt' which would be a linear operator. It should
still work (you would use this instead of AZ_gmres), but I haven't looked at it in years. I'm not sure it is
in any parameter lists or in Stratimikos  ... but I suspect it could be added easily if deemed useful.

-Ray


On 09/26/13 14:04, Erik Boman wrote:
> I don't think one iteration of GMRES is a linear operator because the
> step length is nonlinear. This is true even for CG. It may still work as
> a preconditioner in most cases, but the safer approach would be to use
> FGMRES in the outer iteration.
>
> Erik
>
> Bartlett, Roscoe A. wrote:
>>
>> *Bui,*
>>
>> * *
>>
>> *Let me try to throw an example together which should not be hard if
>> this indeed will work (the only using would be if a single GMRES
>> iteration is not really a linear operator).  If it get it working I
>> will commit it to  the set of Stratimikos examples in Trilinos and
>> then send out the git patch (which should work on the public clone of
>> Trilinos).*
>>
>> * *
>>
>> *-Ross*
>>
>> * *
>>
>>
>>
>> *From:* Hoang Giang Bui [mailto:hgbk2008 at gmail.com]
>> *Sent:* Wednesday, September 25, 2013 7:16 PM
>> *To:* Pawlowski, Roger P; Bartlett, Roscoe A.
>> *Cc:* Heroux, Michael A; Cyr, Eric C; trilinos-users at software.sandia.gov
>> *Subject:* Re: [Trilinos-Users] ilut setting
>>
>>
>>
>>
>> Dear Ross, Pawlowski et al
>>
>> The idea of using NOX is quite simple but is not familiar to me since
>> I didn't use NOX. Additionally, I would want to use AztecOO as
>> sub-preconditioner in block preconditioning. Therefore, I prefer the
>> approach of Ross and wrap it under Teko. However, I didn't quite get
>> the idea since so many interfaces involved. I wrote a snippet to see
>> if I understand this concept correctly:
>>
>>
>> //create the Stratimikos linear solver builder
>> Stratimikos::DefaultLinearSolverBuilder strat;
>> strat.setParameterList(xtraParam); //xtraParam is where I define GMRES
>> and related params
>>
>> //first create a DefaultInverseLinearOp which contains an AztecOO solver
>> Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double>  >
>> PrecOpFactory = strat.createLinearSolveStrategy("AztecOO");
>>
>> Teuchos::RCP<Thyra::LinearOpWithSolveBase<double>  >  PrecOp =
>> PrecOpFactory->createOp();
>>
>> Teuchos::RCP<Thyra::DefaultInverseLinearOp<double>  >  InverseOp =
>> Teuchos::rcp(new Thyra::DefaultInverseLinearOp<double>(PrecOp));
>>
>> //second create a LinearOpWithSolveBase to hold the solver as
>> preconditioner
>> Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double>  >
>> SolveOpFactory = strat.createLinearSolveStrategy("AztecOO"); //I reuse
>> AztecOO here for simplicity
>>
>> Teuchos::RCP<Thyra::LinearOpWithSolveBase<double>  >  SolveOp =
>> SolveOpFactory->createOp();
>>
>> SolveOpFactory->initializePreconditionedOp(....,
>>
>> &*SolveOp,Thyra::SUPPORT_SOLVE_FORWARD_ONLY);
>>
>> I wonder what should be put under dot in initializePreconditionedOp ?
>>
>> Bui
>>
>>
>>
>>
>>
>> On Tue, Sep 24, 2013 at 2:59 PM, Pawlowski, Roger P
>> <rppawlo at sandia.gov<mailto:rppawlo at sandia.gov>>  wrote:
>>
>> Just to follow up on Ross's comment.  The NOX epetra group does allow
>> you to use the aztecoo preconditioners by performing a single GMRES
>> step as he said.  You can get access to this by calling the
>> NOX::Epetra::Group::applyRightPreconditioning() method and selecting
>> the aztecoo preconditioner in the parameter list for building the
>> aztecoo linear system.
>>
>> Roger
>>
>> ------------------------------------------------------------------------
>>
>> *From:* trilinos-users-bounces at software.sandia.gov
>> <mailto:trilinos-users-bounces at software.sandia.gov>
>> [trilinos-users-bounces at software.sandia.gov
>> <mailto:trilinos-users-bounces at software.sandia.gov>] on behalf of
>> Heroux, Michael A [maherou at sandia.gov<mailto:maherou at sandia.gov>]
>> *Sent:* Monday, September 23, 2013 12:13 PM
>> *To:* Bartlett, Roscoe A.; Cyr, Eric C; Hoang Giang Bui
>> *Cc:* trilinos-users at software.sandia.gov
>> <mailto:trilinos-users at software.sandia.gov>
>> *Subject:* Re: [Trilinos-Users] ilut setting
>>
>> Ross' approach is a reasonable path.  AztecOO does have a
>> SetPreconditioner, ConstructPrecondition, and DestroyPreconditioner
>> method.  But it does not have an ApplyPreconditioner.  It would be
>> very straightforward to create such a method using
>> ConstructPreconditioner as the starting point, but I hesitate add code
>> to AztecOO  if the suggested approach that Ross outlined works for you.
>>
>>
>>
>> Mike
>>
>>
>>
>> *From: *<Bartlett>, Ross Bartlett<bartlettra at ornl.gov
>> <mailto:bartlettra at ornl.gov>>
>> *Date: *Monday, September 23, 2013 10:55 AM
>>
>> *To: *Eric Cyr<hgbk2008 at gmail.com<mailto:hgbk2008 at gmail.com>>
>>
>>
>> *Cc: *"trilinos-users at software.sandia.gov
>> <mailto:trilinos-users at software.sandia.gov>>
>> *Subject: *Re: [Trilinos-Users] [EXTERNAL] Re: ilut setting
>>
>>
>>
>> I think NOX uses AztecOO as a precoditioner by doing a single GMRES
>> iteration.  Roger Pawlowski would know more.
>>
>>
>>
>> Stratimikos does not directly support such a usage but you could use
>> Stratimikos and Thyra to construct such a thing.  You would use one
>> Stratikmikos linear solver builder to create an AztecOO
>> LinearOpWtihSolveBase object with GMRES and a single iteration then
>> you would wrap this within a Thyra::DefaultInverseLinearOp object (see
>> http://trilinos.sandia.gov/packages/docs/dev/packages/thyra/doc/html/classThyra_1_1DefaultInverseLinearOp.html
>> ).  You would then use another Stratimikos linear solver build to
>> create a Thyra::LinearOpWithSolveFactoryBase object (i.e. Belos) and
>> then create a new LinearOpWithSolveBase object where you set the
>> AztecOO GMRESS solver set that as a preconditioner using
>> Thyra::LinearOpWithSolveFactoryBase::initializePreconditionedOp() (see
>> http://trilinos.sandia.gov/packages/docs/dev/packages/thyra/doc/html/classThyra_1_1LinearOpWithSolveFactoryBase.html#a3c0e4b3df20383b09eb0f4b671576783
>> ).
>>
>>
>>
>> This is only a few lines of code once you have the two
>> Stratimikos::DefaultLinearSolverBuilder objects created.
>>
>>
>>
>> -Ross
>>
>>
>>
>> *From:* trilinos-users-bounces at software.sandia.gov
>> <mailto:trilinos-users-bounces at software.sandia.gov>
>> [mailto:trilinos-users-bounces at software.sandia.gov] *On Behalf Of
>> *Cyr, Eric C
>> *Sent:* Monday, September 23, 2013 11:02 AM
>> *To:* Hoang Giang Bui
>> *Cc:* trilinos-users at software.sandia.gov
>> <mailto:trilinos-users at software.sandia.gov>
>> *Subject:* Re: [Trilinos-Users] [EXTERNAL] Re: ilut setting
>>
>>
>>
>> Bui
>>
>>
>>
>> So to clarify you would like to use only the preconditioners within
>> AztecOO? Furthermore you want to apply this repeatedly with out
>> rebuilding at each iteration?
>>
>>
>>
>> I'm not sure if only the preconditioners from AztecOO can be used, can
>> someone else answer this?
>>
>>
>>
>> Eric
>>
>>
>> On Sep 23, 2013, at 5:28 AM, "Hoang Giang Bui"<hgbk2008 at gmail.com
>> <mailto:hgbk2008 at gmail.com>>  wrote:
>>
>>      I think the problem lies in line 1044 of
>>      Thyra_AztecOOLinearOpWithSolveFactory.cpp. The
>>      ConstructPreconditioner try to set AZ_pre_calc to AZ_reuse so the
>>      method failed as the previous message showed. However, I'm not
>>      sure that uncomment the next line will do the thing correctly.
>>      Because I tried to uncomment and it can pass the AZ_pre_calc error
>>      but the preconditioner take so long to compute. I think the
>>      Stratimikos's AztecOO interface try to run several iteration steps
>>      to compute the approximation instead of applying the
>>      preconditioner directly. This is not what I want, I think there
>>      should be option to apply AztecOO's preconditioner directly
>>      instead of trial solving steps. This will lead to recomputing the
>>      preconditioner at every iteration steps of outer loop.
>>
>>      Bui
>>
>>      This is the parameter list I used for my test:
>>      <ParameterList name="AztecOO">
>>          <Parameter name="Type" type="string" value="AztecOO"/>
>>          <ParameterList name="Forward Solve">
>>              <ParameterList name="AztecOO Settings">
>>                <Parameter name="Aztec Solver" type="string"
>>      value="BiCGStab"/>
>>                <Parameter name="Aztec Preconditioner" type="string"
>>      value="ilut"/>
>>                <Parameter name="Overlap" type="int" value="1"/>
>>                <Parameter name="Graph Fill" type="int" value="0"/>
>>                <Parameter name="Drop Tolerance" type="double" value="0.0"/>
>>                <Parameter name="Fill Factor" type="double" value="10.0"/>
>>                <Parameter name="Steps" type="int" value="3"/>
>>                <Parameter name="Polynomial Order" type="int" value="3"/>
>>                <Parameter name="RCM Reordering" type="string"
>>      value="Enabled"/>
>>                <Parameter name="Orthogonalization" type="string"
>>      value="Classical"/>
>>                <Parameter name="Size of Krylov Subspace" type="int"
>>      value="1000"/>
>>                <Parameter name="Convergence Test" type="string"
>>      value="r0"/>
>>                <Parameter name="Ill-Conditioning Threshold"
>>      type="double" value="1e+11"/>
>>                <Parameter name="Output Frequency" type="int" value="1"/>
>>              </ParameterList>
>>              <Parameter name="Max Iterations" type="int" value="1"/>
>>              <Parameter name="Tolerance" type="double" value="1e-9"/>
>>          </ParameterList>
>>          <Parameter name="Output Every RHS" type="bool" value="false"/>
>>          <ParameterList name="VerboseObject">
>>              <Parameter name="Output File" type="string" value="none"/>
>>              <Parameter name="Verbosity Level" type="string" value="none"/>
>>          </ParameterList>
>>        </ParameterList>
>>
>>
>>      On 09/20/2013 06:04 PM, Cyr, Eric C wrote:
>>
>>          I don't know why. I've never used the native Aztec
>>          preconditioners and always use IFPACK (though not the ILUT).
>>          Its possible (likely even) that the default AztecOO parameters
>>          specified in Stratimikos don't satisfy the requirements of
>>          Aztec as stated in the error text you sent below. I'd look at
>>          AZ_pre_calc and try to see if its set to "reuse". If so I
>>          think thats the problem.
>>
>>
>>
>>          Eric
>>
>>
>>
>>          *From: *Hoang Giang Bui<hgbk2008 at gmail.com
>>          <mailto:hgbk2008 at gmail.com>>
>>          *Date: *Friday, September 20, 2013 3:04 AM
>>          *To: *Eric Cyr<eccyr at sandia.gov<mailto:eccyr at sandia.gov>>
>>          *Subject: *Re: [EXTERNAL] Re: [Trilinos-Users] ilut setting
>>
>>
>>
>>
>>              Hi again,
>>
>>              I realized that this was not correct since
>>              Thyra::LinearOpBase and Epetra_Operator are not
>>              incompatible. I fixed that using Thyra::block1x1 and it
>>              seems working. However with this setting
>>
>>              <ParameterList name="AztecOO">
>>                  <Parameter name="Type" type="string" value="AztecOO"/>
>>                  <ParameterList name="Forward Solve">
>>                      <ParameterList name="AztecOO Settings">
>>                        <Parameter name="Aztec Solver" type="string"
>>              value="BiCGStab"/>
>>                        <Parameter name="Aztec Preconditioner"
>>              type="string" value="ilut"/>
>>                        <Parameter name="Overlap" type="int" value="1"/>
>>                        <Parameter name="Graph Fill" type="int" value="0"/>
>>                        <Parameter name="Drop Tolerance" type="double"
>>              value="0.0"/>
>>                        <Parameter name="Fill Factor" type="double"
>>              value="10.0"/>
>>                        <Parameter name="Steps" type="int" value="3"/>
>>                        <Parameter name="Polynomial Order" type="int"
>>              value="3"/>
>>                        <Parameter name="RCM Reordering" type="string"
>>              value="Enabled"/>
>>                        <Parameter name="Orthogonalization"
>>              type="string" value="Classical"/>
>>                        <Parameter name="Size of Krylov Subspace"
>>              type="int" value="1000"/>
>>                        <Parameter name="Convergence Test" type="string"
>>              value="r0"/>
>>                        <Parameter name="Ill-Conditioning Threshold"
>>              type="double" value="1e+11"/>
>>                        <Parameter name="Output Frequency" type="int"
>>              value="1"/>
>>                      </ParameterList>
>>                      <Parameter name="Max Iterations" type="int"
>>              value="1"/>
>>                      <Parameter name="Tolerance" type="double"
>>              value="1e-9"/>
>>                  </ParameterList>
>>                  <Parameter name="Output Every RHS" type="bool"
>>              value="true"/>
>>                  <ParameterList name="VerboseObject">
>>                      <Parameter name="Output File" type="string"
>>              value="none"/>
>>                      <Parameter name="Verbosity Level" type="string"
>>              value="default"/>
>>                  </ParameterList>
>>                </ParameterList>
>>
>>
>>
>>              generates a runtime error:
>>
>>
>>              *******************************************************
>>                      ***** Problem: Epetra::CrsMatrix
>>                      ***** Preconditioned BICGSTAB solution
>>                      *****
>>              Thyra::DefaultInverseLinearOp<double>{lows=Thyra::AztecOOLinearOpWithSolve{fwdOp=Thyra::EpetraLinearOp{op='Epetra_CrsMatrix',rangeDim=41443,domainDim=41443}},fwdSolveCriteria=DEFAULT,adjSolveCriteria=DEFAULT}
>>                      ***** No scaling
>>
>>              *******************************************************
>>
>>                              iter:    0           residual = 1.000000e+00
>>
>>
>>              *******************************************************
>>                       ***** Problem: Epetra::CrsMatrix
>>                       ***** Preconditioned BICGSTAB solution
>>                       ***** ILUT( fill-in = 1.000e+01, drop = 0.000e+00)
>>                        ***** with overlap = 1
>>                       ***** No scaling
>>                       ***** NOTE: convergence VARIES when the total
>>              number of
>>                       *****       processors is changed.
>>
>>              *******************************************************
>>
>>               Error:    Did not find previous factorization (requested
>>                   by setting options[AZ_pre_calc] to AZ_reuse).
>>                   To find this factorization, the following
>>                   parameters must match the previous factorization:
>>                        1) Total number of unknowns.
>>                        2) options[AZ_overlap]
>>                        3) options[AZ_scaling]
>>                        4) options[AZ_precond]
>>                        5) options[AZ_reorder]
>>                        6) options[AZ_type_overlap]
>>                        7) options[AZ_subdomain_solve]
>>                        8) options[AZ_graph_fill]
>>                        9) params[AZ_ilut_fill]
>>                       10) params[AZ_drop]
>>                       11) data_org[AZ_name]
>>
>>              Do you know why ?
>>
>>              Bui
>>
>>
>>
>>
>>              On 09/20/2013 12:19 AM, Hoang Giang Bui wrote:
>>
>>
>>                  Hi Eric
>>
>>                  I manage to know all the AztecOO settings by reading
>>                  Stratimikos. However, it still didn't work as I
>>                  expected. Even with the setting "Aztec Solver" as "LU"
>>                  doesn't work as well. That's kind of strange. I want
>>                  to ask if I did program correctly or not. In my case
>>                  for the first step there are no pressure, so matrix
>>                  has just one block. To build preconditioner using
>>                  Teko, I invoke:
>>
>>                  Teuchos::RCP<const Teko::InverseFactory>  inverseA =
>>                  mInverseLibrary->getInverseFactory("AztecOO");
>>
>>                  Teko::LinearOp opA =
>>                  Thyra::epetraLinearOp(blkExtractor.GetBlock(0, 0));
>>
>>                  Teko::LinearOp inv_rA = Teko::buildInverse(*inverseA,
>>                  opA); //build the preconditioner using Teko InverseFactory
>>
>>                  Teuchos::RCP<const Epetra_Operator>  opPrec =
>>                  Teuchos::rcp_dynamic_cast<const Epetra_Operator>(inv_rA);
>>
>>                  Teuchos::RCP<Epetra_Operator>  mPrec =
>>                  Teuchos::rcp_const_cast<Epetra_Operator>(opPrec);
>>
>>                  .............
>>
>>                  aztec_solver.SetPrecOperator(*mPrec);
>>
>>                  Do you think this is correct?
>>
>>                  Bui
>>
>>
>>
>>
>>
>>                  On Thu, Sep 19, 2013 at 8:14 PM, Cyr, Eric C
>>                  <eccyr at sandia.gov<mailto:eccyr at sandia.gov>>  wrote:
>>
>>                  Yeah, thats what you want. Instead of using the string
>>                  "Ifpack" when you
>>                  call the sub solve, use "AztecOO". Teko is "smart"
>>                  enough to know that a
>>                  solver can be used as a preconditioner. (You can also
>>                  use Amesos, if you
>>                  want a direct solver) If you'd like to know how to
>>                  specify a
>>                  preconditioner for Azteco when used as a sub solve let
>>                  me know.
>>
>>                  Eric
>>
>>
>>                  On 9/19/13 12:04 PM, "Hoang Giang Bui"
>>                  <hgbk2008 at gmail.com<mailto:hgbk2008 at gmail.com>>  wrote:
>>
>>                  >Hi Eric
>>                  >
>>                  >Yes, that what I currently use. For example in
>>                  attached file. (The
>>                  >AztecOO section is irrelevant)
>>                  >
>>                  >Bui
>>                  >
>>                  >On 09/19/2013 07:45 PM, Cyr, Eric C wrote:
>>                  >>  Are you using the Teko parameter list, and defining
>>                  an InverseFactory
>>                  >>  through it?
>>                  >>
>>                  >>  Eric
>>                  >>
>>                  >>  On 9/19/13 11:08 AM, "Hoang Giang Bui"
>>                  <hgbk2008 at gmail.com<mailto:hgbk2008 at gmail.com>>  wrote:
>>                  >>
>>                  >>>  Dear Eric
>>                  >>>
>>                  >>>  I think Stratimikos support wrapping AztecOO as
>>                  solver but not
>>                  >>>  preconditioner. Do you have a working script that
>>                  allows AztecOO as
>>                  >>>  preconditioner ?
>>                  >>>
>>                  >>>  Bui
>>                  >>>
>>                  >>>
>>                  >>>  On 09/19/2013 03:42 PM, Cyr, Eric C wrote:
>>                  >>>>  It is already an inverse factory in Teko. AztecOO
>>                  is wrapped by
>>                  >>>>  Stratimikos. Anything in Stratimikos is
>>                  accessible to Teko (This is
>>                  >>>>how
>>                  >>>>  you access Ifpack for instance). If you go to the
>>                  Stratimikos page and
>>                  >>>>  look up the AztecOO parameters you should be able
>>                  to modify those to
>>                  >>>>be
>>                  >>>>  accessible by Teko. Let me know if you need more
>>                  help.
>>                  >>>>
>>                  >>>>  Regarding the ILUT in Ifpack. I have to say that
>>                  I don't use it and
>>                  >>>>  think
>>                  >>>>  there may be problems with its implementation.
>>                  The ILU in Ifpack is
>>                  >>>>  solid
>>                  >>>>  and we get a lot of use out it.
>>                  >>>>
>>                  >>>>  Eric
>>                  >>>>
>>                  >>>>  On 9/19/13 4:43 AM, "Hoang Giang Bui"
>>                  <hgbk2008 at gmail.com<mailto:hgbk2008 at gmail.com>>  wrote:
>>                  >>>>
>>                  >>>>>  That's fine. I can use AztecOO with ilut since
>>                  it works well in my
>>                  >>>>>  case.
>>                  >>>>>  In fact, I used Ifpack as sub-preconditioner for
>>                  my block
>>                  >>>>>  preconditioner. This block preconditioner
>>                  doesn't work since
>>                  >>>>>  Ifpack_ILUT
>>                  >>>>>  failed. It is not sure that AztecOO ILUT will be
>>                  success if it's
>>                  >>>>>  applied
>>                  >>>>>  as sub-preconditioners though. However I wonder
>>                  if there is a way to
>>                  >>>>>  introduce AztecOO as an InverseFactory in Teko?
>>                  >>>>>
>>                  >>>>>  Best regards
>>                  >>>>>  Giang Bui
>>                  >>>>>
>>                  >>>>>
>>                  >>>>>
>>                  >>>>>  On 09/18/2013 07:16 PM, Erik Boman wrote:
>>                  >>>>>>  Unfortunately, the incomplete factorization
>>                  parameters are not well
>>                  >>>>>>  documented so looking at the source code is
>>                  best. Yes, Ifpack_ILUT
>>                  >>>>>>  interprets "ilut level-of-fill" to be the ratio
>>                  of nonzeros in the
>>                  >>>>>>  preconditioner relative to the original matrix,
>>                  so a value of 1.0
>>                  >>>>>>or a
>>                  >>>>>>  bit higher is reasonable. (Note: 1.0 in Ifpack
>>                  roughly corresponds
>>                  >>>>>>to
>>                  >>>>>>  0 in AztecOO.) If you set it to 50, you will in
>>                  most cases compute a
>>                  >>>>>>  complete factorization (but in a very slow
>>                  way!). As a sanity check,
>>                  >>>>>>  you can try using Amesos (by default, KLU) as
>>                  your subdomain solver.
>>                  >>>>>>
>>                  >>>>>>  There are several other parameters in Ifpack
>>                  that may be helpful.
>>                  >>>>>>For
>>                  >>>>>>  example, the absolute and relative threshold
>>                  parameters may make
>>                  >>>>>>your
>>                  >>>>>>  preconditioner more robust. Confusingly, these
>>                  prescribe a diagonal
>>                  >>>>>>  perturbation (aka Manteuffel shift).
>>                  >>>>>>
>>                  >>>>>>  Erik
>>                  >>>>>>
>>                  >>>>>>  Hoang Giang Bui wrote:
>>                  >>>>>>>  Hi Chris
>>                  >>>>>>>
>>                  >>>>>>>  I didn't see the range of value of ilut level
>>                  of fill in Ifpack
>>                  >>>>>>>user
>>                  >>>>>>>  guide yet. But looking in the code of
>>                  Ifpack_ILUT, I'm quite sure
>>                  >>>>>>>  that the level of fill indicates the number of
>>                  nonzeros in the
>>                  >>>>>>>  preconditioner compare to initial nonzeros.
>>                  Therefore its value can
>>                  >>>>>>>  be bigger than 1.
>>                  >>>>>>>
>>                  >>>>>>>  Giang Bui
>>                  >>>>>>>
>>                  >>>>>>>
>>                  >>>>>>>  On 09/18/2013 04:43 AM, Chris Jackson wrote:
>>                  >>>>>>>>  I don¹t know if this applies to the one you
>>                  are talking about,
>>                  >>>>>>>>but I
>>                  >>>>>>>>  recall that with one of the ILUT
>>                  preconditioners the fill value is
>>                  >>>>>>>>  expected to be between 0.0 and 1.0.  1.0
>>                  meant essentially we were
>>                  >>>>>>>>  finding the inverse of the matrix.  Before we
>>                  figured that out we
>>                  >>>>>>>>  were always astounded that BiCGStab was
>>                  converging in 2
>>                  >>>>>>>>iterationsŠ
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>  Chris Jackson
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>  *From:*
>>                  trilinos-users-bounces at software.sandia.gov
>>                  <mailto:trilinos-users-bounces at software.sandia.gov>
>>                  >>>>>>>>
>>                  [mailto:trilinos-users-bounces at software.sandia.gov
>>                  <mailto:trilinos-users-bounces at software.sandia.gov>]
>>                  *On Behalf Of
>>                  >>>>>>>>  *Hoang Giang Bui
>>                  >>>>>>>>  *Sent:* Tuesday, September 17, 2013 6:28 PM
>>                  >>>>>>>>  *To:* Erik Boman
>>                  >>>>>>>>  *Cc:* trilinos-users at software.sandia.gov
>>                  <mailto:trilinos-users at software.sandia.gov>
>>                  >>>>>>>>  *Subject:* Re: [Trilinos-Users] ilut setting
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>  Hi Erik
>>                  >>>>>>>>
>>                  >>>>>>>>  I tried to increase "fact: ilut
>>                  level-of-fill" to 50.0 and it
>>                  >>>>>>>>still
>>                  >>>>>>>>  doesn't converge. In fact, with this value
>>                  the preconditioner take
>>                  >>>>>>>>  very long to compute. Is there other option
>>                  we can set for Ifpack
>>                  >>>>>>>>  ILUT?
>>                  >>>>>>>>
>>                  >>>>>>>>  Giang Bui
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>  On Tue, Sep 17, 2013 at 8:01 PM, Erik Boman
>>                  <egboman at sandia.gov<mailto:egboman at sandia.gov>
>>                  >>>>>>>>  <mailto:egboman at sandia.gov
>>                  <mailto:egboman at sandia.gov>>>  wrote:
>>                  >>>>>>>>
>>                  >>>>>>>>  Giang Bui,
>>                  >>>>>>>>
>>                  >>>>>>>>  You will not get exactly the same answer
>>                  because the ILUT
>>                  >>>>>>>>  implementations are different in AztecOO and
>>                  Ifpack. Also, the
>>                  >>>>>>>>  meaning of "fill" may be slightly different.
>>                  Try increase the
>>                  >>>>>>>>value
>>                  >>>>>>>>  of  "fact: ilut level-of-fill" in the Ifpack
>>                  parameter to see if
>>                  >>>>>>>>you
>>                  >>>>>>>>  can make it converge.
>>                  >>>>>>>>
>>                  >>>>>>>>  Erik
>>                  >>>>>>>>
>>                  >>>>>>>>  Hoang Giang Bui wrote:
>>                  >>>>>>>>
>>                  >>>>>>>>  Dear Trilinos developers,
>>                  >>>>>>>>
>>                  >>>>>>>>  Is there difference between 2 settings:
>>                  >>>>>>>>
>>                  >>>>>>>>  solver_parameters.set("AZ_solver", "AZ_bicgstab")
>>                  >>>>>>>>  solver_parameters.set("AZ_kspace", 1000)
>>                  >>>>>>>>  solver_parameters.set("AZ_output", 100)
>>                  >>>>>>>>  solver_parameters.set("AZ_precond",
>>                  "AZ_dom_decomp")
>>                  >>>>>>>>  solver_parameters.set("AZ_subdomain_solve",
>>                  "AZ_ilut")
>>                  >>>>>>>>  solver_parameters.set("AZ_drop", 0.0)
>>                  >>>>>>>>  solver_parameters.set("AZ_ilut_fill", 10.0)
>>                  >>>>>>>>  aztec_solver.SetParameters(solver_parameter,
>>                  true);
>>                  >>>>>>>>
>>                  >>>>>>>>  And:
>>                  >>>>>>>>  solver_parameters.set("AZ_solver", "AZ_bicgstab")
>>                  >>>>>>>>  solver_parameters.set("AZ_kspace", 1000)
>>                  >>>>>>>>  solver_parameters.set("AZ_output", 100)
>>                  >>>>>>>>  solver_parameters.set("AZ_precond", "AZ_none")
>>                  >>>>>>>>  aztec_solver.SetParameters(solver_parameter,
>>                  true);
>>                  >>>>>>>>  IFPreconditionerType = "ILUT"
>>                  >>>>>>>>  preconditioner_parameters.set("fact: ilut
>>                  level-of-fill", 10.0);
>>                  >>>>>>>>  preconditioner_parameters.set("fact: drop
>>                  tolerance", 0.0);
>>                  >>>>>>>>  overlap_level = 1
>>                  >>>>>>>>  mPrec =
>>                  Teuchos::rcp(Factory.Create(IFPreconditionerType,&rA,
>>                  >>>>>>>>  overlap_level));
>>                  >>>>>>>>  mPrec->SetParameters(preconditioner_parameters);
>>                  >>>>>>>>
>>                  aztec_solver.SetPrecOperator(&*(mpPrec->GetOperator()))
>>                  >>>>>>>>
>>                  >>>>>>>>  I thought it all resulted in
>>                  Additive-Schwartz preconditioner with
>>                  >>>>>>>>  ILUT as local subdomain preconditioner. But
>>                  it gives different
>>                  >>>>>>>>  performance. The first setting works well but
>>                  the second failed to
>>                  >>>>>>>>  converge.
>>                  >>>>>>>>
>>                  >>>>>>>>  Best regards
>>                  >>>>>>>>  Giang Bui
>>                  >>>>>>>>
>>                  >>>>>>>>  _______________________________________________
>>                  >>>>>>>>  Trilinos-Users mailing list
>>                  >>>>>>>>  Trilinos-Users at software.sandia.gov
>>                  <mailto:Trilinos-Users at software.sandia.gov>
>>                  >>>>>>>>  <mailto:Trilinos-Users at software.sandia.gov
>>                  <mailto:Trilinos-Users at software.sandia.gov>>
>>                  >>>>>>>>
>>                  http://software.sandia.gov/mailman/listinfo/trilinos-users
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>
>>                  >>>>>>>>  --
>>                  >>>>>>>>  With Best Regards !
>>                  >>>>>>>>  Giang Bui
>>                  >>>>>>>>
>>                  >>>>>>>>  To learn and to excel
>>                  >>>>>>>>
>>                  >>>>>  _______________________________________________
>>                  >>>>>  Trilinos-Users mailing list
>>                  >>>>>  Trilinos-Users at software.sandia.gov
>>                  <mailto:Trilinos-Users at software.sandia.gov>
>>                  >>>>>
>>                  http://software.sandia.gov/mailman/listinfo/trilinos-users
>>                  >
>>
>>
>>
>>
>>                  --
>>                  With Best Regards !
>>                  Giang Bui
>>
>>                  To learn and to excel
>>
>>
>>
>>
>>
>>
>>
>>
>> --
>> With Best Regards !
>> Giang Bui
>>
>> To learn and to excel
>>
>> ------------------------------------------------------------------------
>>
>> _______________________________________________
>> Trilinos-Users mailing list
>> Trilinos-Users at software.sandia.gov
>> http://software.sandia.gov/mailman/listinfo/trilinos-users
>>
>
> _______________________________________________
> Trilinos-Users mailing list
> Trilinos-Users at software.sandia.gov
> http://software.sandia.gov/mailman/listinfo/trilinos-users

-- 
Ray Tuminaro                      phone: (925) 294-2564
MS 9159                           fax:   (925) 294-2234
Sandia National Laboratories      email: rstumin at sandia.gov
PO Box 969                        http://www.cs.sandia.gov/~rstumin
Livermore, CA 94551




More information about the Trilinos-Users mailing list