<div dir="ltr"><div dir="ltr">On Wed, Jul 14, 2021 at 11:01 AM Tang, Qi <<a href="mailto:tangqi@msu.edu">tangqi@msu.edu</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 style="overflow-wrap: break-word;">
<div>Thanks a lot for the explanation, Matt and Stefano. That helps a lot.</div>
<div><br>
</div>
<div>Just to confirm, the comment in src/ts/impls/implicit/theta/theta.c seems to indicates TS solves U_{n+1}  in its SNES/KSP solve, but it actually solves the update dU_n in U_{n+1} = U_n - lambda*dU_n in the solve. Right?</div>
<div><br>
</div>
<div>It actually makes a lot sense, because KSPSolve in TSSolve reports it uses zero initial guess. So if what I said is true, that effectively means it uses U0 as the initial guess.<br></div></div></blockquote><div><br></div><div>Yes, this makes it fit the SNES API.</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 style="overflow-wrap: break-word;"><div>
Qi</div>
<div><br>
<blockquote type="cite">
<div>On Jul 14, 2021, at 2:56 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div>
<br>
<div>
<div dir="ltr" style="font-family:Helvetica;font-size:16px;font-style:normal;font-variant-caps:normal;font-weight:normal;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none">
<div dir="ltr">On Wed, Jul 14, 2021 at 4:43 AM Stefano Zampini <<a href="mailto:stefano.zampini@gmail.com" target="_blank">stefano.zampini@gmail.com</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>Qi
<div><br>
</div>
<div>Backward Euler is a special case of Theta methods in PETSc (Theta=1). In src/ts/impls/implicit/theta/theta.c on top of SNESTSFormFunction_Theta you have some explanation of what is solved for at each time step (see below). SNES then solves for
 the Newton update dy_n  and the next Newton iterate is computed as x_{n+1} = x_{n} - lambda * dy_n. Hope this helps.</div>
</div>
</blockquote>
<div><br>
</div>
<div>In other words, you should be able to match the initial residual to</div>
<div><br>
</div>
<div>  F(t + dt, 0, -Un / dt)</div>
<div><br>
</div>
<div>for your IFunction. However, it is really not normal to use U = 0. The default is to use U = U0</div>
<div>as the initial guess I think.</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>
<div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">/*</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">  This defines the nonlinear equation that is to be solved with SNES</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">  G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;min-height:24px">
<span style="font-variant-ligatures:no-common-ligatures"></span><br>
</div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">  Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">  otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">  U = (U_{n+1} + U0)/2</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo;color:rgb(64,11,217)">
<span style="font-variant-ligatures:no-common-ligatures">*/</span></div>
<div style="margin:0px;font-stretch:normal;font-size:21px;line-height:normal;font-family:Menlo">
<span style="font-variant-ligatures:no-common-ligatures;color:rgb(47,180,29)">static</span><span style="font-variant-ligatures:no-common-ligatures"><span> </span>PetscErrorCode SNESTSFormFunction_Theta(SNES
 snes,Vec x,Vec y,TS ts)</span></div>
<div><span style="font-variant-ligatures:no-common-ligatures"><br>
</span></div>
<div><br>
<blockquote type="cite">
<div>On Jul 14, 2021, at 6:12 AM, Tang, Qi <<a href="mailto:tangqi@msu.edu" target="_blank">tangqi@msu.edu</a>> wrote:</div>
<br>
<div>
<div>
<div>Hi,</div>
<div><br>
</div>
<div>During the process to experiment the suggestion Matt made, we ran into some questions regarding to TSSolve vs KSPSolve. We got different initial unpreconditioned residual using two solvers. Let’s say we solve the problem with backward Euler and
 there is no rhs. We guess TSSolve solves</div>
<div>(U^{n+1}-U^n)/dt = A U^{n+1}.</div>
<div>(We only provides IJacobian in this case and turn on TS_LINEAR.)</div>
<div>So we guess the initial unpreconditioned residual would be ||U^n/dt||_2, which seems different from the residual we got from a backward Euler stepping we implemented by ourself through KSPSolve.</div>
<div><br>
</div>
<div>Do we have some misunderstanding on TSSolve? </div>
<div><br>
</div>
<div>Thanks,</div>
<div>
<div>Qi</div>
<div>T5@LANL</div>
<div><br>
</div>
<br>
</div>
<div><br>
<blockquote type="cite">
<div>On Jul 7, 2021, at 3:54 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div>
<br>
<div>
<div dir="ltr" style="font-family:Helvetica;font-size:16px;font-style:normal;font-variant-caps:normal;font-weight:normal;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none">
<div dir="ltr">On Wed, Jul 7, 2021 at 2:33 PM Jorti, Zakariae <<a href="mailto:zjorti@lanl.gov" target="_blank">zjorti@lanl.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>
<div id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552divtagdefaultwrapper" dir="ltr" style="font-size:12pt;font-family:Calibri,Helvetica,sans-serif">
<p>Hi Matt,</p>
<p><br>
</p>
<p>Thanks for your quick reply. </p>
<p>I have not completely understood your suggestion, could you please elaborate a bit more? </p>
<p>For your convenience, here is how I am proceeding for the moment in my code: </p>
<p><br>
</p>
<p>TSGetKSP(ts,&ksp);<br>
</p>
<p>KSPGetPC(ksp,&pc);  <br>
</p>
<p>PCSetType(pc,PCFIELDSPLIT);</p>
<p>PCFieldSplitSetDetectSaddlePoint(pc,PETSC_TRUE);</p>
<p>PCSetUp(pc);</p>
<p>PCFieldSplitGetSubKSP(pc, &n, &subksp);</p>
<p>KSPGetPC(subksp[1], &(subpc[1]));</p>
</div>
</div>
</blockquote>
<div>I do not like the two lines above. We should not have to do this. <br>
</div>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
<div>
<div id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552divtagdefaultwrapper" dir="ltr" style="font-size:12pt;font-family:Calibri,Helvetica,sans-serif">
<p><span style="font-size:12pt">KSPSetOperators(subksp[1],T,T);</span></p>
</div>
</div>
</blockquote>
<div> In the above line, I want you to use a separate preconditioning matrix M, instead of T. That way, it will provide<br>
</div>
<div>the preconditioning matrix for your Schur complement problem.</div>
<div><br>
</div>
<div>  Thanks,</div>
<div><br>
</div>
<div>      Matt</div>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
<div>
<div id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552divtagdefaultwrapper" dir="ltr" style="font-size:12pt;font-family:Calibri,Helvetica,sans-serif">
<p><span style="font-size:12pt">KSPSetUp(subksp[1]);</span><br>
</p>
<p>PetscFree(subksp);</p>
<p>TSSolve(ts,X);<br>
</p>
<p><br>
</p>
<p>Thank you.</p>
<p>Best,</p>
<p><br>
</p>
<p>Zakariae</p>
</div>
<hr style="display:inline-block;width:1015.22px">
<div id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552divRplyFwdMsg" dir="ltr">
<font face="Calibri, sans-serif" style="font-size:11pt"><b>From:</b><span> </span>Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
<b>Sent:</b><span> </span>Wednesday, July 7, 2021 12:11:10 PM<br>
<b>To:</b><span> </span>Jorti, Zakariae<br>
<b>Cc:</b><span> </span><a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>; Tang, Qi; Tang, Xianzhu<br>
<b>Subject:</b><span> </span>[EXTERNAL] Re: [petsc-users] Problem with PCFIELDSPLIT</font>
<div> </div>
</div>
<div>
<div dir="ltr">
<div dir="ltr">On Wed, Jul 7, 2021 at 1:51 PM Jorti, Zakariae via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">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 dir="ltr">
<div id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552gmail-m_4230588803823514004divtagdefaultwrapper" dir="ltr" style="font-size:12pt;font-family:Calibri,Helvetica,sans-serif,Helvetica,EmojiFont,"Apple Color Emoji","Segoe UI Emoji",NotoColorEmoji,"Segoe UI Symbol","Android Emoji",EmojiSymbols">
<p>Hi,</p>
<p><br>
</p>
<p>I am trying to build a PCFIELDSPLIT preconditioner for a matrix </p>
<p>J =  [A00  A01]</p>
<p>       [A10  A11] </p>
<p>that has the following shape: </p>
<p><br>
</p>
<p>M_{user}^{-1} = [I   -ksp(A00) A01] [ksp(A00)           0] [I                        0]</p>
<p>                          [0                        I]  [0               ksp(T)] [-A10 ksp(A00)  I ]<br>
</p>
<p><br>
</p>
<p>where T is a user-defined Schur complement approximation that replaces the true Schur complement S:= A11 - A10 ksp(A00) A01.</p>
<p><br>
</p>
<p>I am trying to do something similar to this example (lines 41--45 and 116--121): <a href="https://urldefense.com/v3/__https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex70.c.html__;!!HXCxUKc!hoEfgnaraTfQoSgAiplsc6GJ_HuPXN88m5AJVy1gb7WVMNkGENDnJ3zToOGlhw$" id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552gmail-m_4230588803823514004LPlnk826214" target="_blank">https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex70.c.html</a></p>
<p><br>
</p>
<p>The problem I have is that I manage to replace S with T on a separate single linear system but not for the linear systems generated by my time-dependent PDE. Even if I set the preconditioner <span style="font-size:12pt">M_{user}^{-1}
 correctly, the T matrix gets replaced by S in the preconditioner once I call TSSolve. </span></p>
<p><span style="font-size:12pt">Do you have any suggestions how to fix this knowing that the matrix J does not change over time?</span></p>
<div><span style="font-size:12pt"></span><br>
</div>
</div>
</div>
</blockquote>
<div>I don't like how it is done in that example for this very reason.</div>
<div><br>
</div>
<div>When I want to use a custom preconditioning matrix for the Schur complement, I always give a preconditioning matrix M to the outer solve.</div>
<div>Then PCFIELDSPLIT automatically pulls the correct block from M, (1,1) for the Schur complement, for that preconditioning matrix without</div>
<div>extra code. Can you do this?</div>
<div><br>
</div>
<div>  Thanks,</div>
<div><br>
</div>
<div>    Matt</div>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
<div dir="ltr">
<div id="gmail-m_-7168509508827232524gmail-m_-7613167163141428994gmail-m_4217968377908507552gmail-m_4230588803823514004divtagdefaultwrapper" dir="ltr" style="font-size:12pt;font-family:Calibri,Helvetica,sans-serif,Helvetica,EmojiFont,"Apple Color Emoji","Segoe UI Emoji",NotoColorEmoji,"Segoe UI Symbol","Android Emoji",EmojiSymbols">
<p><span style="font-size:12pt">Many thanks.</span></p>
<p><span style="font-size:12pt"><br>
</span></p>
<p><span style="font-size:12pt">Best regards,</span></p>
<p><span style="font-size:12pt"><br>
</span></p>
<p><span style="font-size:12pt">Zakariae   </span></p>
<div><span style="font-size:12pt"></span><br>
</div>
</div>
</div>
</blockquote>
</div>
<br clear="all">
<div><br>
</div>
--<span> </span><br>
<div dir="ltr">
<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="https://urldefense.com/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!HXCxUKc!hoEfgnaraTfQoSgAiplsc6GJ_HuPXN88m5AJVy1gb7WVMNkGENDnJ3wB3dcMFw$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<br clear="all">
<div><br>
</div>
--<span> </span><br>
<div dir="ltr">
<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="https://urldefense.com/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!HXCxUKc!hoEfgnaraTfQoSgAiplsc6GJ_HuPXN88m5AJVy1gb7WVMNkGENDnJ3wB3dcMFw$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a></div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<br>
</div>
</div>
</blockquote>
</div>
<br>
</div>
</div>
</blockquote>
</div>
<br clear="all">
<div><br>
</div>
--<span> </span><br>
<div dir="ltr">
<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="https://urldefense.com/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!HXCxUKc!msQrz7__TrpOmaTVhvY1yLAlDQXNJ5jcYVAxF4lcpyLrZqt2lFe22bkbuJMizA$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a></div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<br>
</div>

</blockquote></div><br clear="all"><div><br></div>-- <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>