<div dir="ltr"><div dir="ltr">On Wed, Sep 29, 2021 at 6:24 AM Karthikeyan Chockalingam - STFC UKRI <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk">karthikeyan.chockalingam@stfc.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 lang="EN-GB" style="overflow-wrap: break-word;">
<div class="gmail-m_391186430279648065WordSection1">
<p class="MsoNormal"><span>Thank you Mathew. Now, it is all making sense to me.<u></u><u></u></span></p>
<p class="MsoNormal"><span><u></u> <u></u></span></p>
<p class="MsoNormal"><span>From data file </span>
ksp_ex45_N511_gpu_2.txt <u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">KSPSolve (53%) + KSPSetup (0%) = PCSetup (16%) + PCApply (37%).<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">However, you said “So an iteration would mostly consist of MatMult + PCApply, with some vector work”</p></div></div></blockquote><div><br></div><div>1) You do one solve, but 2 KSPSetUp()s. You must be running on more than one process and using Block-Jacobi . Half the time is spent in the solve (53%)</div><div><br></div><div><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px">KSPSetUp               2 1.0 5.3149e-01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 2.4e+01  0  0  0  0  1   0  0  0  0  1     0       0      0 0.00e+00    0 0.00e+00  0
KSPSolve               1 1.0 1.5837e+02 1.1 8.63e+11 1.0 6.8e+02 2.1e+06 4.4e+03 53100100100 95  53100100100 96 10881   11730   1022 6.40e+03 1021 8.17e-03 100
</pre><div><span style="color:rgb(0,0,0);font-family:"Courier New",Courier,monospace,arial,sans-serif;font-size:14px;white-space:pre-wrap"><br></span></div>2) The preconditioner look like BJacobi-ILU. The setup time is 16%, which is all setup of the individual blocks, and this is all used by the numerical ILU factorization.</div><div><br><div><span style="color:rgb(0,0,0);font-family:"Courier New",Courier,monospace,arial,sans-serif;font-size:14px;white-space:pre-wrap">PCSetUp                2 1.0 4.9623e+01 1.3 1.45e+09 1.0 0.0e+00 0.0e+00 0.0e+00 16  0  0  0  0  16  0  0  0  0    58       0      2 6.93e+03    0 0.00e+00  0
PCSetUpOnBlocks        1 1.0 4.9274e+01 1.3 1.45e+09 1.0 0.0e+00 0.0e+00 0.0e+00 15  0  0  0  0  15  0  0  0  0    59       0      2 6.93e+03    0 0.00e+00  0
</span><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px">MatLUFactorNum         1 1.0 4.6126e+01 1.3 1.45e+09 1.0 0.0e+00 0.0e+00 0.0e+00 14  0  0  0  0  14  0  0  0  0    63       0      2 6.93e+03    0 0.00e+00  0
MatILUFactorSym        1 1.0 2.5110e+00 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  1  0  0  0  0   1  0  0  0  0     0       0      0 0.00e+00    0 0.00e+00  0
</pre><div><br></div>3) The preconditioner application takes 37% of the time, which is all solving the factors and recorded in MatSolve(). Matrix multiplication takes 4%.</div><div><br></div><div><span style="color:rgb(0,0,0);font-family:"Courier New",Courier,monospace,arial,sans-serif;font-size:14px;white-space:pre-wrap">PCApply              341 1.0 1.3068e+02 1.6 2.96e+11 1.0 0.0e+00 0.0e+00 0.0e+00 37 34  0  0  0  37 34  0  0  0  4516    4523      1 5.34e+02    0 0.00e+00 100</span><br class="gmail-Apple-interchange-newline"></div><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px">MatSolve             341 1.0 1.3009e+02 1.6 2.96e+11 1.0 0.0e+00 0.0e+00 0.0e+00 36 34  0  0  0  36 34  0  0  0  4536    4538      1 5.34e+02    0 0.00e+00 100</pre><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px">MatMult              341 1.0 1.0774e+01 1.1 2.96e+11 1.0 6.9e+02 2.1e+06 2.0e+00  4 34100100  0   4 34100100  0 54801   66441      2 5.86e+03    0 0.00e+00 100
<br></pre>4) The significant vector time is all in norms (11%) since they are really slow on the GPU.</div><div><br><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px">VecNorm              342 1.0 6.2261e+01129.9 4.57e+10 1.0 0.0e+00 0.0e+00 6.8e+02 11  5  0  0 15  11  5  0  0 15  1466   196884      0 0.00e+00  342 2.74e-03 100</pre><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px">VecTDot              680 1.0 1.7107e+00 1.3 9.09e+10 1.0 0.0e+00 0.0e+00 1.4e+03  1 10  0  0 29   1 10  0  0 29 106079   133922      0 0.00e+00  680 5.44e-03 100
VecAXPY              681 1.0 3.2036e+00 1.7 9.10e+10 1.0 0.0e+00 0.0e+00 0.0e+00  1 11  0  0  0   1 11  0  0  0 56728   58367    682 5.34e+02    0 0.00e+00 100
VecAYPX              339 1.0 2.6502e+00 1.8 4.53e+10 1.0 0.0e+00 0.0e+00 0.0e+00  1  5  0  0  0   1  5  0  0  0 34136   34153    339 2.71e-03    0 0.00e+00 100</pre><pre class="gmail-aLF-aPX-K0-aPE" style="font-family:"Courier New",Courier,monospace,arial,sans-serif;margin-top:0px;margin-bottom:0px;white-space:pre-wrap;color:rgb(0,0,0);font-size:14px"><br></pre>So the solve time is:</div><div><br></div><div>  53% ~ 37% + 4% + 11%</div><div><br></div><div>and the setup time is about 16%. I was wrong about the SetUp time being included, as it is outside the event:</div><div><br></div><div>  <a href="https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/ksp/interface/itfunc.c#L852">https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/ksp/interface/itfunc.c#L852</a></div><div><br></div><div>It looks like the remainder of the time (23%) is spent preallocating the matrix.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div><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 lang="EN-GB" style="overflow-wrap: break-word;"><div class="gmail-m_391186430279648065WordSection1">
<p class="MsoNormal">The MalMult event is 4 %. How does this event figure into the above equation; if preconditioning (MatMult + PCApply) is included in KSPSolve?<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Best,<u></u><u></u></p>
<p class="MsoNormal">Karthik.<u></u><u></u></p>
<p class="MsoNormal"><span><u></u> <u></u></span></p>
<div style="border-right:none;border-bottom:none;border-left:none;border-top:1pt solid rgb(181,196,223);padding:3pt 0cm 0cm">
<p class="MsoNormal"><b><span style="font-size:12pt;color:black">From: </span></b><span style="font-size:12pt;color:black">Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
<b>Date: </b>Wednesday, 29 September 2021 at 10:58<br>
<b>To: </b>"Chockalingam, Karthikeyan (STFC,DL,HC)" <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>><br>
<b>Cc: </b>Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>>, "<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>" <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>><br>
<b>Subject: </b>Re: [petsc-users] %T (percent time in this phase)<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<div>
<p class="MsoNormal">On Wed, Sep 29, 2021 at 5:52 AM Karthikeyan Chockalingam - STFC UKRI <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>> wrote:<u></u><u></u></p>
</div>
<div>
<blockquote style="border-top:none;border-right:none;border-bottom:none;border-left:1pt solid rgb(204,204,204);padding:0cm 0cm 0cm 6pt;margin-left:4.8pt;margin-right:0cm">
<div>
<div>
<p class="MsoNormal">That was helpful. I would like to provide some additional details of my run on cpus and gpus. Please find the following attachments:<u></u><u></u></p>
<p class="MsoNormal"> <u></u><u></u></p>
<ol start="1" type="1">
<li class="gmail-m_391186430279648065gmail-m-3409338445930824232msolistparagraph">
graph.pdf a plot showing overall time and various petsc events.<u></u><u></u></li><li class="gmail-m_391186430279648065gmail-m-3409338445930824232msolistparagraph">
ksp_ex45_N511_cpu_6.txt data file of the log_summary <u></u><u></u></li><li class="gmail-m_391186430279648065gmail-m-3409338445930824232msolistparagraph">
ksp_ex45_N511_gpu_2.txt data file of the log_summary <u></u><u></u></li></ol>
<p class="MsoNormal"> <u></u><u></u></p>
<p class="MsoNormal">I used the following petsc options for cpu<u></u><u></u></p>
<p class="MsoNormal"><span style="color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black">mpirun -n 6 ./ex45 -log_summary -da_grid_x 511 -da_grid_y 511 -da_grid_z 511 -dm_mat_type mpiaij -dm_vec_type mpi -ksp_type cg -pc_type bjacobi -ksp_monitor</span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black">and for gpus</span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black">mpirun -n 1 ./ex45 -log_summary -da_grid_x 511 -da_grid_y 511 -da_grid_z 511  -dm_mat_type mpiaijcusparse -dm_vec_type mpicuda -ksp_type cg -pc_type bjacobi
 -ksp_monitor</span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black">to run the following problem</span><u></u><u></u></p>
<p class="MsoNormal"><span style="color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black"><a href="https://petsc.org/release/src/ksp/ksp/tutorials/ex45.c.html" target="_blank">https://petsc.org/release/src/ksp/ksp/tutorials/ex45.c.html</a></span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black">From the above code, I see is there no individual function called
</span>KSPSetUp(), so I gather KSPSetDM, KSPSetComputeInitialGuess, KSPSetComputeRHS, kSPSetComputeOperators all are timed together as KSPSetUp.<span style="font-size:12pt;color:black"> For this example, is KSPSetUp time and KSPSolve time mutually exclusive?</span><u></u><u></u></p>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">No, KSPSetUp() will be contained in KSPSolve() if it is called automatically.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<blockquote style="border-top:none;border-right:none;border-bottom:none;border-left:1pt solid rgb(204,204,204);padding:0cm 0cm 0cm 6pt;margin-left:4.8pt;margin-right:0cm">
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt;color:black">In your response you said that</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal">   “PCSetUp() time may be in KSPSetUp() or it maybe in PCApply() it depends on how much of the preconditioner construction can take place early, so depends exactly on the preconditioner
 used.”<u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black">I don’t see a explicit call to PCSetUp() or  PCApply() in ex45; so for this particular preconditioner (bjacobi) how can I tell how they
 are timed?</span><u></u><u></u></p>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">They are all inside KSPSolve(). If you have a preconditioned linear solve, the oreconditioning happens during the iteration. So an iteration would mostly<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">consist of MatMult + PCApply, with some vector work.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<blockquote style="border-top:none;border-right:none;border-bottom:none;border-left:1pt solid rgb(204,204,204);padding:0cm 0cm 0cm 6pt;margin-left:4.8pt;margin-right:0cm">
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt;color:black">I am hoping to time KSP solving and preconditioning mutually exclusively.</span><u></u><u></u></p>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">I am not sure that concept makes sense here. See above.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">  Thanks,<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">     Matt<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<blockquote style="border-top:none;border-right:none;border-bottom:none;border-left:1pt solid rgb(204,204,204);padding:0cm 0cm 0cm 6pt;margin-left:4.8pt;margin-right:0cm">
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black">Kind regards,</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black">Karthik.</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:12pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"> <u></u><u></u></p>
<div style="border-right:none;border-bottom:none;border-left:none;border-top:1pt solid rgb(181,196,223);padding:3pt 0cm 0cm">
<p class="MsoNormal"><b><span style="font-size:12pt;color:black">From:
</span></b><span style="font-size:12pt;color:black">Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>><br>
<b>Date: </b>Tuesday, 28 September 2021 at 19:19<br>
<b>To: </b>"Chockalingam, Karthikeyan (STFC,DL,HC)" <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>><br>
<b>Cc: </b>"<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>" <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>><br>
<b>Subject: </b>Re: [petsc-users] %T (percent time in this phase)</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<p class="MsoNormal"> <u></u><u></u></p>
<div>
<p class="MsoNormal" style="margin-bottom:12pt"><u></u> <u></u></p>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal">On Sep 28, 2021, at 12:11 PM, Karthikeyan Chockalingam - STFC UKRI <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>>
 wrote:<u></u><u></u></p>
</div>
<p class="MsoNormal"> <u></u><u></u></p>
<div>
<div>
<p class="MsoNormal">Thanks for Barry for your response.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">I was just benchmarking the problem with various preconditioner on cpu and gpu. I understand, it is not possible to get mutually exclusive timing.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">However, can you tell if KSPSolve time includes both PCSetup and PCApply? And if KSPSolve and KSPSetup are mutually exclusive? Likewise for PCSetUp and PCApply.<u></u><u></u></p>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<p class="MsoNormal">   If you do not call KSPSetUp() separately from KSPSolve() then its time is included with KSPSolve().<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">   PCSetUp() time may be in KSPSetUp() or it maybe in PCApply() it depends on how much of the preconditioner construction can take place early, so depends exactly on the preconditioner
 used.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">   So yes the answer is not totally satisfying. The one thing I would recommend is to not call KSPSetUp() directly and then KSPSolve() will always include the total time of the
 solve plus all setup time. PCApply will contain all the time to apply the preconditioner but may also include some setup time.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  Barry<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">Best,<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">Karthik.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div style="border-right:none;border-bottom:none;border-left:none;border-top:1pt solid rgb(181,196,223);padding:3pt 0cm 0cm">
<div>
<p class="MsoNormal"><b><span style="font-size:12pt">From:<span class="gmail-m_391186430279648065gmail-m-3409338445930824232apple-converted-space"> </span></span></b><span style="font-size:12pt">Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>><br>
<b>Date:<span class="gmail-m_391186430279648065gmail-m-3409338445930824232apple-converted-space"> </span></b>Tuesday, 28 September 2021 at 16:56<br>
<b>To:<span class="gmail-m_391186430279648065gmail-m-3409338445930824232apple-converted-space"> </span></b>"Chockalingam, Karthikeyan (STFC,DL,HC)" <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>><br>
<b>Cc:<span class="gmail-m_391186430279648065gmail-m-3409338445930824232apple-converted-space"> </span></b>"<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>" <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>><br>
<b>Subject:<span class="gmail-m_391186430279648065gmail-m-3409338445930824232apple-converted-space"> </span></b>Re: [petsc-users] %T (percent time in this phase)</span><u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<div>
<p class="MsoNormal" style="margin-bottom:12pt"><br>
<br>
<u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<div>
<p class="MsoNormal">On Sep 28, 2021, at 10:55 AM, Karthikeyan Chockalingam - STFC UKRI <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>>
 wrote:<u></u><u></u></p>
</div>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt">Hello,</span><u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt"> </span><u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt">I ran ex45 in the KPS tutorial, which is a 3D finite-difference Poisson problem. I noticed from the output from using the flag -log_summary that for
 various events their respective %T (percent time in this phase) do not add up to 100 but rather exceeds 100. So, I gather there is some overlap among these events. I am primarily looking at the events KSPSetUp, KSPSolve, PCSetUp and PCSolve. Is it possible
 to get a mutually exclusive %T or Time for these individual events? I have attached  the log_summary output file from my run for your reference.</span><u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt"> </span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<p class="MsoNormal">  For nested solvers it is tricky to get the times to be mutually exclusive because some parts of the building of the preconditioner is for some preconditioners delayed until the
 solve has started. <u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal">  It looks like you are using the default preconditioner options which for this example are taking more or less no time since so many iterations are needed. It is best to use -pc_type
 mg to use geometric multigrid on this problem.<u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal">  Barry<u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal" style="margin-bottom:12pt"><br>
<br>
<u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt">Thanks!</span><u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:12pt">Karthik.</span><u></u><u></u></p>
</div>
</div>
<div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
</div>
<div>
<p class="MsoNormal"><span style="font-size:6pt;font-family:Helvetica">This email and any attachments are intended solely for the use of the named recipients. If you are not the intended recipient
 you must not use, disclose, copy or distribute this email or any of its attachments and should notify the sender immediately and delete this email from your system. UK Research and Innovation (UKRI) has taken every reasonable precaution to minimise risk of
 this email or any attachments containing viruses or malware but the recipient should carry out its own virus and malware checks before opening the attachments. UKRI does not accept any liability for any losses or damages which the recipient may sustain due
 to presence of any viruses. </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><ksp_ex45_N511_cpu_6.txt><u></u><u></u></p>
</div>
</div>
</blockquote>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal" style="margin-bottom:12pt"><u></u> <u></u></p>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><br clear="all">
<u></u><u></u></p>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<p class="MsoNormal">-- <u></u><u></u></p>
<div>
<div>
<div>
<div>
<div>
<div>
<div>
<p class="MsoNormal">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<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><u></u><u></u></p>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</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>