<div dir="ltr"><div dir="ltr"><div dir="ltr">How about I want to determine the ST type on runtime? <div><br></div><div> mpirun -n 1 ./ex3  -eps_type jd -st_ksp_type gmres  -st_pc_type none   -eps_view  -eps_target  0 -eps_monitor  -st_ksp_monitor <br></div><div><br></div><div>ST is indeed STPrecond, but the passed preconditioning matrix is still ignored.</div><div><br></div><div><div>EPS Object: 1 MPI processes</div><div>  type: jd</div><div>    search subspace is orthogonalized</div><div>    block size=1</div><div>    type of the initial subspace: non-Krylov</div><div>    size of the subspace after restarting: 6</div><div>    number of vectors after restarting from the previous iteration: 1</div><div>    threshold for changing the target in the correction equation (fix): 0.01</div><div>  problem type: symmetric eigenvalue problem</div><div>  selected portion of the spectrum: closest to target: 0. (in magnitude)</div><div>  number of eigenvalues (nev): 1</div><div>  number of column vectors (ncv): 17</div><div>  maximum dimension of projected problem (mpd): 17</div><div>  maximum number of iterations: 1700</div><div>  tolerance: 1e-08</div><div>  convergence test: relative to the eigenvalue</div><div>BV Object: 1 MPI processes</div><div>  type: svec</div><div>  17 columns of global length 100</div><div>  vector orthogonalization method: classical Gram-Schmidt</div><div>  orthogonalization refinement: if needed (eta: 0.7071)</div><div>  block orthogonalization method: GS</div><div>  doing matmult as a single matrix-matrix product</div><div>DS Object: 1 MPI processes</div><div>  type: hep</div><div>  solving the problem with: Implicit QR method (_steqr)</div><div>ST Object: 1 MPI processes</div><div>  type: precond</div><div>  shift: 0.</div><div>  number of matrices: 1</div><div>  KSP Object: (st_) 1 MPI processes</div><div>    type: gmres</div><div>      restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>      happy breakdown tolerance 1e-30</div><div>    maximum iterations=90, initial guess is zero</div><div>    tolerances:  relative=0.0001, absolute=1e-50, divergence=10000.</div><div>    left preconditioning</div><div>    using PRECONDITIONED norm type for convergence test</div><div>  PC Object: (st_) 1 MPI processes</div><div>    type: none</div><div>    linear system matrix = precond matrix:</div><div>    Mat Object: 1 MPI processes</div><div><span style="background-color:rgb(255,0,0)">      type: shell</span></div><div>      rows=100, cols=100</div><div> Solution method: jd</div></div><div><br></div><div><br></div><div>Preconding matrix should be a SeqAIJ not shell.</div><div><br></div><div><br></div><div>Fande,</div></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Tue, Nov 5, 2019 at 9:07 AM Jose E. Roman <<a href="mailto:jroman@dsic.upv.es">jroman@dsic.upv.es</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex">Currently, the function that passes the preconditioner matrix is specific of STPRECOND, so you have to add<br>
  ierr = STSetType(st,STPRECOND);CHKERRQ(ierr);<br>
before<br>
  ierr = STPrecondSetMatForPC(st,B);CHKERRQ(ierr);<br>
otherwise this latter call is ignored.<br>
<br>
We may be changing a little bit the way in which ST is initialized, and maybe we modify this as well. It is not decided yet.<br>
<br>
Jose<br>
<br>
<br>
> El 5 nov 2019, a las 0:28, Fande Kong <<a href="mailto:fdkong.jd@gmail.com" target="_blank">fdkong.jd@gmail.com</a>> escribió:<br>
> <br>
> Thanks Jose,<br>
> <br>
> I think I understand now. Another question: what is the right way to setup a linear preconditioning matrix for the inner linear solver of JD?<br>
> <br>
> I was trying to do something like this:<br>
> <br>
>   /*<br>
>      Create eigensolver context<br>
>   */<br>
>   ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);<br>
> <br>
>   /*<br>
>      Set operators. In this case, it is a standard eigenvalue problem<br>
>   */<br>
>   ierr = EPSSetOperators(eps,A,NULL);CHKERRQ(ierr);<br>
>   ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);<br>
>   ierr = EPSGetST(eps,&st);CHKERRQ(ierr);<br>
>   ierr = STPrecondSetMatForPC(st,B);CHKERRQ(ierr);<br>
> <br>
>   /*<br>
>      Set solver parameters at runtime<br>
>   */<br>
>   ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);<br>
> <br>
>   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -<br>
>                       Solve the eigensystem<br>
>      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */<br>
> <br>
>   ierr = EPSSolve(eps);CHKERRQ(ierr);<br>
> <br>
> <br>
> But did not work. A complete example is attached.  I could try to dig into the code, but you may already know the answer.<br>
> <br>
> <br>
> On Wed, Oct 23, 2019 at 3:58 AM Jose E. Roman <<a href="mailto:jroman@dsic.upv.es" target="_blank">jroman@dsic.upv.es</a>> wrote:<br>
> Yes, it is confusing. Here is the explanation: when you use a target, the preconditioner is built from matrix A-sigma*B. By default, instead of TARGET_MAGNITUDE we set LARGEST_MAGNITUDE, and in Jacobi-Davidson we treat this case by setting sigma=PETSC_MAX_REAL. In this case, the preconditioner is built from matrix B. The thing is that in a standard eigenproblem we have B=I, and hence there is no point in using a preconditioner, that is why we set PCNONE.<br>
> <br>
> Jose<br>
> <br>
> <br>
> > El 22 oct 2019, a las 19:57, Fande Kong via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>> escribió:<br>
> > <br>
> > Hi All,<br>
> > <br>
> > It looks like the preconditioner is hard-coded in the Jacobi-Davidson solver. I could not select a preconditioner rather than the default setting.<br>
> > <br>
> > For example, I was trying to select LU, but PC NONE was still used.  I ran standard example 2 in slepc/src/eps/examples/tutorials, and had the following results.<br>
> > <br>
> > <br>
> > Thanks,<br>
> > <br>
> > Fande<br>
> > <br>
> > <br>
> > ./ex2 -eps_type jd -st_ksp_type gmres  -st_pc_type lu   -eps_view  <br>
> > <br>
> > 2-D Laplacian Eigenproblem, N=100 (10x10 grid)<br>
> > <br>
> > EPS Object: 1 MPI processes<br>
> >   type: jd<br>
> >     search subspace is orthogonalized<br>
> >     block size=1<br>
> >     type of the initial subspace: non-Krylov<br>
> >     size of the subspace after restarting: 6<br>
> >     number of vectors after restarting from the previous iteration: 1<br>
> >     threshold for changing the target in the correction equation (fix): 0.01<br>
> >   problem type: symmetric eigenvalue problem<br>
> >   selected portion of the spectrum: largest eigenvalues in magnitude<br>
> >   number of eigenvalues (nev): 1<br>
> >   number of column vectors (ncv): 17<br>
> >   maximum dimension of projected problem (mpd): 17<br>
> >   maximum number of iterations: 1700<br>
> >   tolerance: 1e-08<br>
> >   convergence test: relative to the eigenvalue<br>
> > BV Object: 1 MPI processes<br>
> >   type: svec<br>
> >   17 columns of global length 100<br>
> >   vector orthogonalization method: classical Gram-Schmidt<br>
> >   orthogonalization refinement: if needed (eta: 0.7071)<br>
> >   block orthogonalization method: GS<br>
> >   doing matmult as a single matrix-matrix product<br>
> > DS Object: 1 MPI processes<br>
> >   type: hep<br>
> >   solving the problem with: Implicit QR method (_steqr)<br>
> > ST Object: 1 MPI processes<br>
> >   type: precond<br>
> >   shift: 1.79769e+308<br>
> >   number of matrices: 1<br>
> >   KSP Object: (st_) 1 MPI processes<br>
> >     type: gmres<br>
> >       restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br>
> >       happy breakdown tolerance 1e-30<br>
> >     maximum iterations=90, initial guess is zero<br>
> >     tolerances:  relative=0.0001, absolute=1e-50, divergence=10000.<br>
> >     left preconditioning<br>
> >     using PRECONDITIONED norm type for convergence test<br>
> >   PC Object: (st_) 1 MPI processes<br>
> >     type: none<br>
> >     linear system matrix = precond matrix:<br>
> >     Mat Object: 1 MPI processes<br>
> >       type: shell<br>
> >       rows=100, cols=100<br>
> >  Solution method: jd<br>
> > <br>
> >  Number of requested eigenvalues: 1<br>
> >  Linear eigensolve converged (1 eigenpair) due to CONVERGED_TOL; iterations 20<br>
> >  ---------------------- --------------------<br>
> >             k             ||Ax-kx||/||kx||<br>
> >  ---------------------- --------------------<br>
> >         7.837972            7.71944e-10<br>
> >  ---------------------- --------------------<br>
> > <br>
> > <br>
> > <br>
> <br>
> <ex3.c><br>
<br>
</blockquote></div>