<div dir="ltr"><div dir="ltr">On Thu, Apr 20, 2023 at 6:13 AM Karthikeyan Chockalingam - STFC UKRI 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="msg4059408548263094035">
<div lang="EN-GB" style="overflow-wrap: break-word;">
<div class="m_-6455518359002065503WordSection1">
<p class="MsoNormal">Hello,</p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">I created a new thread, thought would it be more appropriate (and is a continuation of my previous post). I want to construct the below K matrix (which is composed of submatrices)
</p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal"><span style="color:black">K = [A P^T</span><span style="font-size:10pt;color:rgb(33,33,33)"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"> P 0]<u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black">Where K is of type MatMPIAIJ. I first constructed the top left [A] using MatSetValues().<u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black">Now, I would like to construct the bottom left [p] and top right [p^T]
</span>using MatSetValuesLocal().</p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">To use MatSetValuesLocal(), I first have to create a local-to-global mapping using ISLocalToGlobalMappingCreate. I have created two mapping row_mapping and column_mapping.</p></div></div></div></blockquote><div><br></div><div>I do not understand why they are not the same map. Maybe I was unclear before. It looks like you have two fields, say phi and lambda, where lambda is a Lagrange multiplier imposing some constraint. Then you get a saddle point like this. You can imagine matrices</div><div><br></div><div> (phi, phi) --> A</div><div> (phi, lambda) --> P^T</div><div> (lambda, phi) --> P</div><div><br></div><div>So you make a L2G map for the phi field and the lambda field. Oh, you are calling them row and col map, but they are my phi and lambda</div><div>maps. I do not like the row and col names since in P they reverse.</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="msg4059408548263094035"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_-6455518359002065503WordSection1">
<p class="MsoNormal"><span style="color:black">Q1) At what point should I declare
</span>MatSetLocalToGlobalMapping – is it just before I use MatSetValuesLocal()?</p></div></div></div></blockquote><div><br></div><div>Okay, it is good you are asking this because my thinking was somewhat confused. I think the precise steps are:</div><div><br></div><div> 1) Create the large saddle point matrix K</div><div><br></div><div> 1a) We must call <a href="https://petsc.org/main/manualpages/Mat/MatSetLocalToGlobalMapping/">https://petsc.org/main/manualpages/Mat/MatSetLocalToGlobalMapping/</a> on it. In the simplest case, this just maps</div><div> the local rows numbers [0, Nrows) to the global rows numbers [rowStart, rowStart + Nrows).</div><div><br></div><div> 2) To form each piece:</div><div><br></div><div> 2a) Extract that block using <a href="https://petsc.org/main/manualpages/Mat/MatGetLocalSubMatrix/">https://petsc.org/main/manualpages/Mat/MatGetLocalSubMatrix/</a></div><div><br></div><div> This gives back a Mat object that you subsequently restore using <a href="https://petsc.org/main/manualpages/Mat/MatRestoreLocalSubMatrix/">https://petsc.org/main/manualpages/Mat/MatRestoreLocalSubMatrix/</a></div><div><br></div><div> 2b) Insert values using <a href="https://petsc.org/main/manualpages/Mat/MatSetValuesLocal/">https://petsc.org/main/manualpages/Mat/MatSetValuesLocal/</a></div><div><br></div><div> The local indices used for insertion here are indices relative to the block itself, and the L2G map for this matrix</div><div> has been rewritten to insert into that block in the larger matrix. Thus this looks like just calling MatSetValuesLocal()</div><div> on the smaller matrix block, but inserts correctly into the larger matrix.</div><div><br></div><div>Therefore, the code you write code in 2) could work equally well making the large matrix from 1), or independent smaller matrix blocks.</div><div><br></div><div>Does this make sense?</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="msg4059408548263094035"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_-6455518359002065503WordSection1">
<p class="MsoNormal">I will use MatSetLocalToGlobalMapping(K, row_mapping, column_mapping) to build
<span style="color:black">the bottom left [P].<u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black">Q2) Can now I reset the mapping as
</span>MatSetLocalToGlobalMapping(K, column_mapping, row_mapping) to build <span style="color:black">
the top right [P^T]? <u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Many thanks!<span style="color:black"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black">Kind regards,<u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black">Karthik.<u></u><u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="color:black"><u></u> <u></u></span></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal"><span style="font-size:10pt;color:rgb(33,33,33)"><u></u> <u></u></span></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>