<div dir="ltr"><div dir="ltr">On Mon, Aug 21, 2023 at 4:06 AM Jonas Lundgren via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="msg-6373068259469628963">





<div lang="EN-US">
<div class="m_-6373068259469628963WordSection1">
<p class="MsoNormal"><span lang="SV">Dear PETSc users,<u></u><u></u></span></p>
<p class="MsoNormal"><span lang="SV"><u></u> <u></u></span></p>
<p class="MsoNormal">I have a problem regarding the setting of initial guess to KSP when using PCREDISTRIBUTE as the preconditioner. (The reason to why I use PCREDISTRIBUTE is because I have a lot of fixed DOF in my problem, and PCREDISTRIBUTE successfully
 reduces the problem size and therefore speeds up the solving).<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">First, some details:<u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>I use a version of PETSc 3.19.1<u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>The KSP I use is KSPPREONLY, as suggested in the manual pages of PCREDISTRIBUTE:
<a href="https://petsc.org/release/manualpages/PC/PCREDISTRIBUTE/" target="_blank">https://petsc.org/release/manualpages/PC/PCREDISTRIBUTE/</a><u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>I use KSPBCGSL as sub-KSP. I can perfectly well solve my problem using this as my main KSP, but the performance is much worse than when using it as my sub-KSP (under KSPPREONLY+PCREDISTRIBUTE) due to the amount of fixed DOF<u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>I am first solving a state problem, then an adjoint problem using the same linear operator.<u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>The adjoint vector is used as sensitivity information to update a design. After the design update, the state+adjoint problems are solved again with a slightly updated linear operator. This is done for hundreds of (design) iteration steps<u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>I want the initial guess for the state problem to be the state solution from the previous (design) iteration, and same for the adjoint problem<u></u><u></u></p>
<p class="m_-6373068259469628963MsoListParagraph"><u></u><span>-<span style="font:7pt "Times New Roman"">         
</span></span><u></u>I am aware of the default way of setting a custom initial guess: KSPSetInitialGuessNonzero(ksp, PETSC_TRUE) together with providing the actual guess in the x vector in the call to KSPSolve(ksp, b, x)<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">The main problem is that PCREDISTRIBUTE internally doesn’t use the input solution vector (x) when calling KSPSolve() for the sub-KSP. It zeroes out the solution vector (x) when starting to build x = diag(A)^{-1} b in the beginning of PCApply_Redistribute(),
 and uses red->x as the solution vector/initial guess when calling KSPSolve(). Therefore, I cannot reach the sub-KSP with an initial guess.<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Additionally, KSPPREONLY prohibits the use of having a nonzero initial guess (the error message says “it doesn’t make sense”). I guess I can remove the line raising this error and recompile the PETSc libraries, but it still won’t solve
 the previously mentioned problem, which seems to be the hard nut to crack.<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">So far, I have found out that if I create 2 KSP object, one each for the state and adjoint problems, it is enough with calling KSPSetInitialGuessNonzero(subksp, PETSC_TRUE) on the subksp. It seems as if the variable red->x in PCApply_Redistribute()
 is kept untouched in memory between calls to the main KSP and therefore is used as (non-zero) initial guess to the sub-KSP. This has been verified by introducing PetscCall(PetscObjectCompose((PetscObject)pc,"redx",(PetscObject)red->x)); in PCApply_Redistribute(),
 recompiling the PETSc library, and then inserting a corresponding PetscObjectQuery((PetscObject)pc, "redx", (PetscObject *)&redx); in my own program source file.<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">However, I would like to only create 1 KSP to be used with both the state and adjoint problems (same linear operator), for memory reasons. </p></div></div></div></blockquote><div><br></div><div>I  would like to understand this better. Are you trying to save the temporary vector memory for another solver? For BCGSL, this is only a few vectors. It should be much less than the matrices.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="msg-6373068259469628963"><div lang="EN-US"><div class="m_-6373068259469628963WordSection1"><p class="MsoNormal">When I do this, the initial guesses are mixed up between the two problems: the initial guess for
 the adjoint problem is the solution to the state problem in the current design iteration, and the initial guess for the state problem is the solution to the adjoint problem in the previous design iteration. These are very bad guesses that increases the time
 to convergence in each state/adjoint solve.<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">So, the core of the problem (as far as I can understand) is that I want to control the initial guess red->x in PCApply_Redistribute().<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">The only solution I can think of is to include a call to PetscObjectQuery() in PCApply_Redistribute() to obtain a vector with the initial guess from my main program. And then I need to keep track of the initial guesses for my both problems
 in my main program myself (minor problem). This is maybe not the neatest way, and I do not know if this approach affects the performance negatively? Maybe one call each to PetscObjectQuery() and PetscObjectCompose() per call to PCApply_Redistribute() is negligible?
<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Is there another (and maybe simpler) solution to this problem? Maybe I can SCATTER_FORWARD the input x vector in PCApply_Redistribute() before it is zeroed out, together with allowing for non-zero initial guess in KSPPREONLY?<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Any help is welcome!<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Best regards,<u></u><u></u></p>
<p class="MsoNormal">Jonas Lundgren<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
</div>

</div></blockquote></div><br clear="all"><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>