<div dir="ltr"><div dir="ltr">On Tue, Sep 22, 2020 at 12:20 PM Mebratu Wakeni <<a href="mailto:Mebratu.Wakeni@glasgow.ac.uk">Mebratu.Wakeni@glasgow.ac.uk</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 style="font-family:Calibri,Arial,Helvetica,sans-serif;font-size:12pt;color:rgb(0,0,0);background-color:rgb(255,255,255)">
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
Hello,</p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<br>
</p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
I have the following issue regarding TS and field split. </p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<br>
</p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
In a TS setting, I have a linear problem with the IJacobian of the IFunction with a block form<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
| K<span>  <span> </span></span>0 |<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
| 0<span>  <span> </span></span>D |<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
Where K is sparse, while D is a diagonal. The sizes of the square matrices K and D are dynamic (varying between time steps). </p></div></div></blockquote><div><br></div><div>How are you doing this in TS? Are you using just TSStep and calling TSReset() after each step?</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 dir="ltr"><div style="font-family:Calibri,Arial,Helvetica,sans-serif;font-size:12pt;color:rgb(0,0,0);background-color:rgb(255,255,255)"><p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">The fact that, K and D are completely decoupled, and D is just a diagonal motivated us to use fieldsplit techniques rather than a monolithic
 one. However, this needs to recreate a fieldsplit preconditioner at each time step using the indices of K. To achieve this, I use the<span> </span><span style="font-size:10pt;font-family:Menlo">TSSetPreStep(ts, fun) </span>function
 to perform the<span> </span><span style="font-size:10pt;font-family:Menlo">fun(TS ts)</span><span> </span>routine which does following<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<u></u> <u></u></p>
<p style="margin:0cm 0cm 0cm 36pt;font-size:12pt;font-family:Calibri,sans-serif">
<span><span><span>1.<span style="font:7pt "Times New Roman"">    <span> </span></span></span></span></span>Get
 indices of K,<span> </span><u></u> <u></u></p>
<p style="font-family:Calibri,sans-serif;margin:0cm 0cm 0cm 36pt;font-size:medium;color:rgb(0,0,0)">
<u></u> <u></u></p>
<p style="font-family:Calibri,sans-serif;margin:0cm 0cm 0cm 36pt;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo">PetscInt *idx_k; // size_of_k<u></u> <u></u></span></p>
<p style="font-family:Calibri,sans-serif;margin:0cm 0cm 0cm 36pt;font-size:medium;color:rgb(0,0,0)">
<u></u> <u></u></p>
<p style="margin:0cm 0cm 0cm 36pt;font-size:12pt;font-family:Calibri,sans-serif">
<span><span><span>2.<span style="font:7pt "Times New Roman"">    <span> </span></span></span></span></span>Create
 an IS corresponding to K (IS is_k) using<span> </span><u></u> <u></u><br>
<br>
<u></u> <u></u><br>
<br>
<span> </span><span style="font-size:11pt;font-family:Menlo">ISCreatGeneral(mpiComm, size_of_k, idx_k, PETSC_COPY_VALUES, &is_k);<span> </span><u></u> <u></u></span><br>
<br>
<u></u> <u></u></p>
<p style="margin:0cm 0cm 0cm 36pt;font-size:12pt;font-family:Calibri,sans-serif">
<span><span>3.<span style="font:7pt "Times New Roman"">    <span> </span></span></span></span>Then set a fieldsplit preconditioner PC
 pc using is_k through<span> </span><u></u> <u></u><br>
<br>
<u></u> <u></u><br>
<br>
<span style="font-size:11pt;font-family:Menlo">PCFieldSplitSetIS(pc, NULL, is_k);<u></u> <u></u></span><br>
<br>
<span> </span><u></u> <u></u></p>
<p style="margin:0cm 0cm 0cm 36pt;font-size:12pt;font-family:Calibri,sans-serif">
<span><span>4.<span style="font:7pt "Times New Roman"">    <span> </span></span></span></span>Finally the is_k is properly destroyed
 through<span> </span><u></u> <u></u><br>
<br>
<u></u> <u></u></p>
<p style="font-family:Calibri,sans-serif;margin:0cm 0cm 0cm 36pt;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo">ISDestroy(&is_lo);<u></u> </span></p></div></div></blockquote><div><br></div><div>That all looks fine.</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 dir="ltr"><div style="font-family:Calibri,Arial,Helvetica,sans-serif;font-size:12pt;color:rgb(0,0,0);background-color:rgb(255,255,255)">
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
With the above setup, I have also used the following runtime options<span> </span><u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<u></u> <u></u></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ts_monitor<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ts_exact_final_time stepover<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ts_dt 0.1<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ts_max_time 1<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ts_adapt_type none<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black"><u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ksp_type gmres<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-pc_type fieldsplit<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-pc_fieldsplit_diag_use_amat</span></p></div></div></blockquote><div><br></div><div>You do not need this unless you have a separate preconditioner matrix.</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 dir="ltr"><div style="font-family:Calibri,Arial,Helvetica,sans-serif;font-size:12pt;color:rgb(0,0,0);background-color:rgb(255,255,255)"><p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)"><span style="font-size:11pt;font-family:Menlo;color:black"><u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-pc_fieldsplit_type additive<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-fieldsplit_0_ksp_type preonly<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-fieldsplit_0_pc_type lu<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-fieldsplit_0_pc_factor_mat_solver_type mumps<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-fieldsplit_1_ksp_type preonly<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-fieldsplit_1_pc_type jacobi<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black"><u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-snes_atol 1e-8<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-snes_rtol 1e-12<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-ksp_monitor<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-snes_monitor<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="font-size:11pt;font-family:Menlo;color:black">-snes_lag_jacobian 1</span></p></div></div></blockquote><div><br></div><div>You should not need this since your SNES should converge in 1 iterate with an exact linear solver and a linear problem.</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 dir="ltr"><div style="font-family:Calibri,Arial,Helvetica,sans-serif;font-size:12pt;color:rgb(0,0,0);background-color:rgb(255,255,255)">
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="color:black">Result: with the setup described above and the options, I expect to get ksp convergence in a single step (at every time step), however, what I am getting in actuality is that the fieldsplit works as expected only for the first timestep,
 after that it either take too many ksp steps to converge or diverge completely. I have a suspicion that, the way the fieldsplit set up at each time using the TSSetPreStep() somehow adds new split over previous splits. For this, I have also tried using TSSetPostStep()
 function to reset both the ksp and pc as follows<u></u> </span></p></div></div></blockquote><div><br></div><div>First check -snes_view to see if you have the solver you think you do, and also</div><div><br></div><div>  -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -snes_converged_reason</div><div><br></div><div>It sounds to me like something is not being recreated correctly.</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 dir="ltr"><div style="font-family:Calibri,Arial,Helvetica,sans-serif;font-size:12pt;color:rgb(0,0,0);background-color:rgb(255,255,255)">
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="color:black">KSPRest(ksp);<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="color:black">PCReset(pc);<u></u> <u></u></span></p>
<p class="MsoNormal" style="font-family:Calibri,sans-serif;margin:0cm;font-size:medium;color:rgb(0,0,0)">
<span style="color:black"><u></u> <u></u></span></p>
<span style="font-size:12pt;font-family:Calibri,sans-serif;color:black">But these don’t help either.<span> </span><span> </span><span> </span></span><span style="font-family:-webkit-standard;font-size:medium;color:rgb(0,0,0);display:inline"></span><br>
</div>
</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>