<div dir="ltr"><br><div class="gmail_extra"><br><div class="gmail_quote">On Thu, Oct 13, 2016 at 8:21 PM, Barry Smith <span dir="ltr"><<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><br>
Fande,<br>
<br>
What SNES method are you using? If you use SNESKSPONLY I think it is ok, it will solve for the norm minimizing least square solution during the one KSPSolve() and then return.<br></blockquote><div><br></div><div>The problem we are currently working on is a linear problem, but it could be extended to be nonlinear. Yes, you are right. "ksponly" indeed works, and returns the right solution. But the norm of residual still could confuse users because it is not close to zero.<br></div><div><br> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<br>
Yes, if you use SNESNEWTONLS or others though the SNES solver will, as you say, think that progress has not been made.<br>
<br>
I do not like what you propose to do, changing the right hand side of the system the user provides is a nasty and surprising side effect.<br></blockquote><div><br></div><div>I do not like this way either. The reason I posted this code here is that I want to let you know what are inconsistent between the nonlinear solvers and the linear solvers. <br></div><div><br> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<br>
What is your goal? To make it look like the SNES system has had a residual norm reduction?<br></blockquote><div><br></div><div>Yes, I would like to make SNES have a residual reduction. Possibly, we could add something in the converged_test function? For example, the residual vector is temporarily subtracted when evaluating the residual norm if the system has a null space?<br></div><div><br> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<br>
We could generalize you question and ask what about solving for nonlinear problems: find the minimal norm solution of min_x || F(x) - b||. This may or may not belong in Tao, currently SNES doesn't do any kind of nonlinear least squares.<br></blockquote><div><br><br></div><div>It would be great, if we could add this kind of solvers. Tao does have one, I think. I would like to contribute something like this latter (of course, if you are ok with this algorithm), when we are moving to nonlinear problems in our applications. <br><br></div><div>Fande Kong,<br></div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<span class="HOEnZb"><font color="#888888"><br>
Barry<br>
</font></span><div class="HOEnZb"><div class="h5"><br>
<br>
> On Oct 13, 2016, at 5:20 PM, Kong, Fande <<a href="mailto:fande.kong@inl.gov">fande.kong@inl.gov</a>> wrote:<br>
><br>
> One more question.<br>
><br>
> Suppose that we are solving the singular linear system Ax = b. N(A) is the null space of A, and N(A^T) is the null space of the transpose of A.<br>
><br>
> The linear system is solved using SNES, that is, F(x) = Ax-b = Ax -b_r - b_n. Here b_n in N(A^T), and b_r in R(A). During each nonlinear iteration, a linear system A \delta x = F(x) is solved. N(A) is applied to Krylov space during the linear iterating. Before the actual solve "(*ksp->ops->solve)(ksp)" for \delta x, a temporary copy of F(x) is made, F_tmp. N(A^T) is applied to F_tmp. We will get a \delta x. F(x+\delta x ) = A(x+\delta x)-b_r - b_n.<br>
><br>
> F(x+\delta x ) always contain the vector b_n, and then the algorithm never converges because the normal of F is at least 1.<br>
><br>
> Should we apply N(A^T) to F instead of F_tmp so that b_n can be removed from F?<br>
><br>
> MatGetTransposeNullSpace(pmat,<wbr>&nullsp);<br>
> if (nullsp) {<br>
> VecDuplicate(ksp->vec_rhs,&<wbr>btmp);<br>
> VecCopy(ksp->vec_rhs,btmp);<br>
> MatNullSpaceRemove(nullsp,<wbr>btmp);<br>
> vec_rhs = ksp->vec_rhs;<br>
> ksp->vec_rhs = btmp;<br>
> }<br>
><br>
> should be changed to<br>
><br>
> MatGetTransposeNullSpace(pmat,<wbr>&nullsp);<br>
> if (nullsp) {<br>
> MatNullSpaceRemove(nullsp,ksp-<wbr>>vec_rhs);<br>
> }<br>
> ???<br>
><br>
> Or other solutions to this issue?<br>
><br>
><br>
> Fande Kong,<br>
><br>
><br>
><br>
><br>
><br>
> On Thu, Oct 13, 2016 at 8:23 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br>
> On Thu, Oct 13, 2016 at 9:06 AM, Kong, Fande <<a href="mailto:fande.kong@inl.gov">fande.kong@inl.gov</a>> wrote:<br>
><br>
><br>
> On Wed, Oct 12, 2016 at 10:21 PM, Jed Brown <<a href="mailto:jed@jedbrown.org">jed@jedbrown.org</a>> wrote:<br>
> Barry Smith <<a href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a>> writes:<br>
> > I would make that a separate routine that the users would call first.<br>
><br>
> We have VecMDot and VecMAXPY. I would propose adding<br>
><br>
> VecQR(PetscInt nvecs,Vec *vecs,PetscScalar *R);<br>
><br>
> (where R can be NULL).<br>
><br>
> What does R mean here?<br>
><br>
> It means the coefficients of the old basis vectors in the new basis.<br>
><br>
> Matt<br>
><br>
> If nobody working on this, I will be going to take a try.<br>
><br>
> Fande,<br>
><br>
><br>
> Does anyone use the "Vecs" type?<br>
><br>
><br>
><br>
><br>
> --<br>
> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
> -- Norbert Wiener<br>
><br>
<br>
</div></div></blockquote></div><br></div></div>