[petsc-users] (Sub) KSP initial guess with PCREDISTRIBUTE
Matthew Knepley
knepley at gmail.com
Mon Aug 21 10:47:54 CDT 2023
On Mon, Aug 21, 2023 at 4:06 AM Jonas Lundgren via petsc-users <
petsc-users at mcs.anl.gov> wrote:
> Dear PETSc users,
>
>
>
> 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).
>
>
>
> First, some details:
>
> - I use a version of PETSc 3.19.1
>
> - The KSP I use is KSPPREONLY, as suggested in the manual pages
> of PCREDISTRIBUTE:
> https://petsc.org/release/manualpages/PC/PCREDISTRIBUTE/
>
> - 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
>
> - I am first solving a state problem, then an adjoint problem
> using the same linear operator.
>
> - 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
>
> - 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
>
> - 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)
>
>
>
> 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.
>
>
>
> 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.
>
>
>
> 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.
>
>
>
> 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.
>
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.
Thanks,
Matt
> 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.
>
>
>
> 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().
>
>
>
> 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?
>
>
>
> 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?
>
>
>
> Any help is welcome!
>
>
>
>
>
> Best regards,
>
> Jonas Lundgren
>
>
>
--
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230821/210054f2/attachment-0001.html>
More information about the petsc-users
mailing list