<html>
  <head>
    <meta content="text/html; charset=utf-8" http-equiv="Content-Type">
  </head>
  <body bgcolor="#FFFFFF" text="#000000">
    Mark, <br>
      Can you give me the full option line that you want me to use?  I
    currently have:<br>
    <br>
    -ksp_type cg -ksp_monitor -ksp_chebyshev_esteig_random -log_view
    -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1
    -options_left<br>
    <br>
    -sanjay<br>
    <br>
    <div class="moz-cite-prefix">On 5/22/16 2:29 PM, Mark Adams wrote:<br>
    </div>
    <blockquote
cite="mid:CADOhEh48JzX4Av3iSDebE44vBB9JkXMCLo4akGBXs3YFyeMpNQ@mail.gmail.com"
      type="cite">
      <div dir="ltr">Humm, maybe we have version mixup:
        <div><br>
          <div>src/ksp/ksp/impls/cheby/cheby.c:    ierr =
            PetscOptionsBool("-ksp_chebyshev_esteig_random","Use random
            right hand side for
estimate","KSPChebyshevEstEigSetUseRandom",cheb->userandom,&cheb->userandom<br>
            <div><br>
            </div>
            <div>Also, you should use CG. These other options are the
              defaults but CG is not:</div>
            <div><br>
            </div>
            <div>-mg_levels_esteig_ksp_type cg </div>
            <div>-mg_levels_esteig_ksp_max_it 10 </div>
            <div>-mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05<br>
            </div>
          </div>
        </div>
        <div><br>
        </div>
        <div>Anyway. you can also run with -info, which will be very
          noisy, but just grep for GAMG and send me that.</div>
        <div><br>
        </div>
        <div>Mark</div>
        <div><br>
        </div>
        <div><br>
        </div>
      </div>
      <div class="gmail_extra"><br>
        <div class="gmail_quote">On Sat, May 21, 2016 at 6:03 PM, Sanjay
          Govindjee <span dir="ltr"><<a moz-do-not-send="true"
              href="mailto:s_g@berkeley.edu" target="_blank">s_g@berkeley.edu</a>></span>
          wrote:<br>
          <blockquote class="gmail_quote" style="margin:0 0 0
            .8ex;border-left:1px #ccc solid;padding-left:1ex">
            <div bgcolor="#FFFFFF" text="#000000"> Mark,<br>
                I added the option you mentioned but it seems not to use
              it; -options_left reports:<br>
              <br>
              #PETSc Option Table entries:<br>
              -ksp_chebyshev_esteig_random<br>
              -ksp_monitor<br>
              -ksp_type cg<br>
              -log_view<br>
              -options_left<br>
              -pc_gamg_agg_nsmooths 1<br>
              -pc_gamg_type agg<br>
              -pc_type gamg<br>
              #End of PETSc Option Table entries<br>
              There is one unused database option. It is:<br>
              Option left: name:-ksp_chebyshev_esteig_random (no value)
              <div>
                <div class="h5"><br>
                  <br>
                  <br>
                  <div>On 5/21/16 12:36 PM, Mark Adams wrote:<br>
                  </div>
                  <blockquote type="cite">
                    <div dir="ltr">Barry, this is probably the Chebyshev
                      problem.
                      <div><br>
                      </div>
                      <div>Sanjay, this is fixed but has not yet been
                        moved to the master branch.  You can fix this
                        now with with -ksp_chebyshev_esteig_random. 
                        This should recover v3.5 semantics.</div>
                      <div><br>
                      </div>
                      <div>Mark </div>
                      <div class="gmail_extra"><br>
                        <div class="gmail_quote">On Thu, May 19, 2016 at
                          2:42 PM, Barry Smith <span dir="ltr"><<a
                              moz-do-not-send="true"
                              href="mailto:bsmith@mcs.anl.gov"
                              target="_blank"><a class="moz-txt-link-abbreviated" href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a></a>></span>
                          wrote:<br>
                          <blockquote class="gmail_quote"
                            style="margin:0 0 0 .8ex;border-left:1px
                            #ccc solid;padding-left:1ex"><br>
                               We see this occasionally, there is
                            nothing in the definition of GAMG that
                            guarantees a positive definite
                            preconditioner even if the operator was
                            positive definite so we don't think this is
                            a bug in the code. We've found using a
                            slightly stronger smoother, like one more
                            smoothing step seems to remove the problem.<br>
                            <br>
                               Barry<br>
                            <div>
                              <div><br>
                                > On May 19, 2016, at 1:07 PM, Sanjay
                                Govindjee <<a moz-do-not-send="true"
                                  href="mailto:s_g@berkeley.edu"
                                  target="_blank">s_g@berkeley.edu</a>>
                                wrote:<br>
                                ><br>
                                > I am trying to solve a very
                                ordinary nonlinear elasticity problem<br>
                                > using -ksp_type cg -pc_type gamg in
                                PETSc 3.7.0, which worked fine<br>
                                > in PETSc 3.5.3.<br>
                                ><br>
                                > The problem I am seeing is on my
                                first Newton iteration, the Ax=b<br>
                                > solve returns with and Indefinite
                                Preconditioner error
                                (KSPGetConvergedReason == -8):<br>
                                > (log_view.txt output also attached)<br>
                                ><br>
                                >   0 KSP Residual norm
                                8.411630828687e-02<br>
                                >   1 KSP Residual norm
                                2.852209578900e-02<br>
                                >   NO CONVERGENCE REASON: 
                                Indefinite Preconditioner<br>
                                >   NO CONVERGENCE REASON: 
                                Indefinite Preconditioner<br>
                                ><br>
                                > On the next and subsequent Newton
                                iterations, I see perfectly normal<br>
                                > behavior and the problem converges
                                quadratically.  The results look fine.<br>
                                ><br>
                                > I tried the same problem with
                                -pc_type jacobi as well as super-lu, and
                                mumps<br>
                                > and they all work without
                                complaint.<br>
                                ><br>
                                > My run line for GAMG is:<br>
                                > -ksp_type cg -ksp_monitor -log_view
                                -pc_type gamg -pc_gamg_type agg
                                -pc_gamg_agg_nsmooths 1 -options_left<br>
                                ><br>
                                > The code flow looks like:<br>
                                ><br>
                                > ! If no matrix allocation yet<br>
                                > if(Kmat.eq.0) then<br>
                                >   call
                                MatCreate(PETSC_COMM_WORLD,Kmat,ierr)<br>
                                >   call
                                MatSetSizes(Kmat,numpeq,numpeq,PETSC_DETERMINE,PETSC_DETERMINE,ierr)<br>
                                >   call
                                MatSetBlockSize(Kmat,nsbk,ierr)<br>
                                >   call MatSetFromOptions(Kmat,
                                ierr)<br>
                                >   call MatSetType(Kmat,MATAIJ,ierr)<br>
                                >   call
MatMPIAIJSetPreallocation(Kmat,PETSC_NULL_INTEGER,mr(np(246)),PETSC_NULL_INTEGER,mr(np(247)),ierr)<br>
                                >   call
                                MatSeqAIJSetPreallocation(Kmat,PETSC_NULL_INTEGER,mr(np(246)),ierr)<br>
                                > endif<br>
                                ><br>
                                > call MatZeroEntries(Kmat,ierr)<br>
                                ><br>
                                > ! Code to set values in matrix<br>
                                ><br>
                                > call MatAssemblyBegin(Kmat,
                                MAT_FINAL_ASSEMBLY, ierr)<br>
                                > call MatAssemblyEnd(Kmat,
                                MAT_FINAL_ASSEMBLY, ierr)<br>
                                > call
                                MatSetOption(Kmat,MAT_NEW_NONZERO_LOCATIONS,PETSC_TRUE,ierr)<br>
                                ><br>
                                > ! If no rhs allocation yet<br>
                                > if(rhs.eq.0) then<br>
                                >   call VecCreate       
                                (PETSC_COMM_WORLD, rhs, ierr)<br>
                                >   call VecSetSizes      (rhs,
                                numpeq, PETSC_DECIDE, ierr)<br>
                                >   call VecSetFromOptions(rhs, ierr)<br>
                                > endif<br>
                                ><br>
                                > ! Code to set values in RHS<br>
                                ><br>
                                > call VecAssemblyBegin(rhs, ierr)<br>
                                > call VecAssemblyEnd(rhs, ierr)<br>
                                ><br>
                                > if(kspsol_exists) then<br>
                                >   call KSPDestroy(kspsol,ierr)<br>
                                > endif<br>
                                ><br>
                                > call KSPCreate(PETSC_COMM_WORLD,
                                kspsol   ,ierr)<br>
                                > call KSPSetOperators(kspsol, Kmat,
                                Kmat, ierr)<br>
                                > call KSPSetFromOptions(kspsol,ierr)<br>
                                > call KSPGetPC(kspsol, pc ,     
                                ierr)<br>
                                ><br>
                                > call
                                PCSetCoordinates(pc,ndm,numpn,hr(np(43)),ierr)<br>
                                ><br>
                                > call KSPSolve(kspsol, rhs, sol,
                                ierr)<br>
                                > call
                                KSPGetConvergedReason(kspsol,reason,ierr)<br>
                                ><br>
                                > ! update solution, go back to the
                                top<br>
                                ><br>
                                > reason is coming back as -8 on my
                                first Ax=b solve and 2 or 3 after that<br>
                                > (with gamg).  With the other
                                solvers it is coming back as 2 or 3 for<br>
                                > iterative options and 4 if I use
                                one of the direct solvers.<br>
                                ><br>
                                > Any ideas on what is causing the
                                Indefinite PC on the first iteration
                                with GAMG?<br>
                                ><br>
                                > Thanks in advance,<br>
                                > -sanjay<br>
                                ><br>
                                > --<br>
                                >
                                -----------------------------------------------<br>
                                > Sanjay Govindjee, PhD, PE<br>
                                > Professor of Civil Engineering<br>
                                ><br>
                                > 779 Davis Hall<br>
                                > University of California<br>
                                > Berkeley, CA 94720-1710<br>
                                ><br>
                                > Voice:  <a moz-do-not-send="true"
                                  href="tel:%2B1%20510%20642%206060"
                                  value="+15106426060" target="_blank">+1
                                  510 642 6060</a><br>
                                > FAX:    <a moz-do-not-send="true"
                                  href="tel:%2B1%20510%20643%205264"
                                  value="+15106435264" target="_blank">+1
                                  510 643 5264</a><br>
                                ><br>
                                > <a moz-do-not-send="true"
                                  href="mailto:s_g@berkeley.edu"
                                  target="_blank">s_g@berkeley.edu</a><br>
                                > <a moz-do-not-send="true"
                                  href="http://www.ce.berkeley.edu/%7Esanjay"
                                  rel="noreferrer" target="_blank">http://www.ce.berkeley.edu/~sanjay</a><br>
                                ><br>
                                >
                                -----------------------------------------------<br>
                                ><br>
                                > Books:<br>
                                ><br>
                                > Engineering Mechanics of Deformable<br>
                                > Solids: A Presentation with
                                Exercises<br>
                                ><br>
                                > <a moz-do-not-send="true"
href="http://www.oup.com/us/catalog/general/subject/Physics/MaterialsScience/?view=usa&ci=9780199651641"
                                  rel="noreferrer" target="_blank">http://www.oup.com/us/catalog/general/subject/Physics/MaterialsScience/?view=usa&ci=9780199651641</a><br>
                                > <a moz-do-not-send="true"
                                  href="http://ukcatalogue.oup.com/product/9780199651641.do"
                                  rel="noreferrer" target="_blank">http://ukcatalogue.oup.com/product/9780199651641.do</a><br>
                                > <a moz-do-not-send="true"
                                  href="http://amzn.com/0199651647"
                                  rel="noreferrer" target="_blank">http://amzn.com/0199651647</a><br>
                                ><br>
                                ><br>
                                > Engineering Mechanics 3 (Dynamics)
                                2nd Edition<br>
                                ><br>
                                > <a moz-do-not-send="true"
                                  href="http://www.springer.com/978-3-642-53711-0"
                                  rel="noreferrer" target="_blank">http://www.springer.com/978-3-642-53711-0</a><br>
                                > <a moz-do-not-send="true"
                                  href="http://amzn.com/3642537111"
                                  rel="noreferrer" target="_blank">http://amzn.com/3642537111</a><br>
                                ><br>
                                ><br>
                                > Engineering Mechanics 3,
                                Supplementary Problems: Dynamics<br>
                                ><br>
                                > <a moz-do-not-send="true"
                                  href="http://www.amzn.com/B00SOXN8JU"
                                  rel="noreferrer" target="_blank">http://www.amzn.com/B00SOXN8JU</a><br>
                                ><br>
                                ><br>
                                >
                                -----------------------------------------------<br>
                                ><br>
                              </div>
                            </div>
                            > <log_view.txt><br>
                            <br>
                          </blockquote>
                        </div>
                        <br>
                      </div>
                    </div>
                  </blockquote>
                  <br>
                  <pre cols="72">-- 
-----------------------------------------------
Sanjay Govindjee, PhD, PE
Professor of Civil Engineering

779 Davis Hall
University of California
Berkeley, CA 94720-1710

Voice:  <a moz-do-not-send="true" href="tel:%2B1%20510%20642%206060" value="+15106426060" target="_blank">+1 510 642 6060</a>
FAX:    <a moz-do-not-send="true" href="tel:%2B1%20510%20643%205264" value="+15106435264" target="_blank">+1 510 643 5264</a>
<a moz-do-not-send="true" href="mailto:s_g@berkeley.edu" target="_blank">s_g@berkeley.edu</a>
<a moz-do-not-send="true" href="http://www.ce.berkeley.edu/%7Esanjay" target="_blank">http://www.ce.berkeley.edu/~sanjay</a>
-----------------------------------------------

Books:  

Engineering Mechanics of Deformable 
Solids: A Presentation with Exercises
<a moz-do-not-send="true" href="http://www.oup.com/us/catalog/general/subject/Physics/MaterialsScience/?view=usa&ci=9780199651641" target="_blank">http://www.oup.com/us/catalog/general/subject/Physics/MaterialsScience/?view=usa&ci=9780199651641</a>
<a moz-do-not-send="true" href="http://ukcatalogue.oup.com/product/9780199651641.do" target="_blank">http://ukcatalogue.oup.com/product/9780199651641.do</a>
<a moz-do-not-send="true" href="http://amzn.com/0199651647" target="_blank">http://amzn.com/0199651647</a>

Engineering Mechanics 3 (Dynamics) 2nd Edition
<a moz-do-not-send="true" href="http://www.springer.com/978-3-642-53711-0" target="_blank">http://www.springer.com/978-3-642-53711-0</a>
<a moz-do-not-send="true" href="http://amzn.com/3642537111" target="_blank">http://amzn.com/3642537111</a>

Engineering Mechanics 3, Supplementary Problems: Dynamics 
<a moz-do-not-send="true" href="http://www.amzn.com/B00SOXN8JU" target="_blank">http://www.amzn.com/B00SOXN8JU</a>

-----------------------------------------------
</pre>
                </div>
              </div>
            </div>
          </blockquote>
        </div>
        <br>
      </div>
    </blockquote>
    <br>
  </body>
</html>