<div dir="ltr"><div dir="ltr">On Mon, Apr 24, 2023 at 10:22 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 class="msg-968874379040273580">
<div lang="EN-GB" style="overflow-wrap: break-word;">
<div class="m_-968874379040273580WordSection1">
<p class="MsoNormal"><span style="font-size:11pt">Hello,<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">I was able to construct the below K matrix (using submatrices P and P^T), which is of type MATAIJ
<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black">K = [A P^T</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> P 0]<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">and solved them using a direct solver.</span></p></div></div></div></blockquote><div><br></div><div>I modified your example to create either AIJ or Nest matrices, but use the same assembly code:</div><div><br></div><div> <a href="https://gitlab.com/petsc/petsc/-/merge_requests/6368">https://gitlab.com/petsc/petsc/-/merge_requests/6368</a></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="msg-968874379040273580"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_-968874379040273580WordSection1"><p class="MsoNormal"><span style="font-size:11pt">
<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">However, I was reading online that this is a saddle point problem and I should be employing PCFIELDSPLIT.
<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">Since I have one monolithic matrix K, I was not sure how to split the fields.</span></p></div></div></div></blockquote><div><br></div><div>With this particular matrix, you can use</div><div><br></div><div> -pc_fieldsplit_detect_saddle_point</div><div><br></div><div>and it will split it automatically.</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="msg-968874379040273580"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_-968874379040273580WordSection1">
<p class="MsoNormal"><span style="font-size:11pt">Best regards,<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">Karthik.<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt"><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" style="margin-bottom:12pt"><b><span style="font-size:12pt;color:black">From:
</span></b><span style="font-size:12pt;color:black">Chockalingam, Karthikeyan (STFC,DL,HC) <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>><br>
<b>Date: </b>Wednesday, 19 April 2023 at 17:52<br>
<b>To: </b>Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</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] Setting up a matrix for Lagrange multiplier<u></u><u></u></span></p>
</div>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> I have declared the mapping
</span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> ISLocalToGlobalMapping mapping;</span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, n, nindices, PETSC_COPY_VALUES, &mapping);</span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black">But when I use MatSetValuesLocal(), how do I know the above mapping is employed because it is not one of the parameters passed to the function?</span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black">Thank you.</span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black">Kind regards,</span><u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:11pt;color:black">Karthik.</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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" style="margin-bottom:12pt"><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>Tuesday, 18 April 2023 at 16:21<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] Setting up a matrix for Lagrange multiplier</span><u></u><u></u></p>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">On Tue, Apr 18, 2023 at 11:16 AM Karthikeyan Chockalingam - STFC UKRI <</span><a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank"><span style="font-size:11pt">karthikeyan.chockalingam@stfc.ac.uk</span></a><span style="font-size:11pt">>
wrote:</span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="m_-968874379040273580m-4128250693350498318p1">Thank you for your response. I spend some time understanding how<u></u><u></u></p>
<p class="m_-968874379040273580m-4128250693350498318p1">MatSetValuesLocal and ISLocalToGlobalMappingCreate work.<u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You can look at SNES ex28 where we do this with DMCOMPOSITE.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">Q1) Will the matrix K be of type MATMPIAIJ or MATIS?</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black">K = [A P^T</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> P 0]</span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">I assume MPIAIJ since IS is only used for Neumann-Neumann decompositions.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">Q2) Can I use both MatSetValues() to MatSetValuesLocal() to populate K? Since I have already used MatSetValues() to construct A.</span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You can, and there would be no changes in serial if K is exactly the upper left block, but in parallel global indices would change.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">Q3) What are the advantages of using MatSetValuesLocal()? Is it that I can construct P directly using local indies and map the entrees to the global
index in K?</span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You have a monolithic K, so that you can use sparse direct solvers to check things. THis is impossible with separate storage.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">Q4) I probably don’t have to construct an independent P matrix</span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You wouldn't in this case.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Thanks,</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Matt</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">Best regards,</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt">Karthik.</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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" style="margin-bottom:12pt"><b><span style="font-size:12pt;color:black">From:
</span></b><span style="font-size:12pt;color:black">Matthew Knepley <</span><a href="mailto:knepley@gmail.com" target="_blank"><span style="font-size:12pt">knepley@gmail.com</span></a><span style="font-size:12pt;color:black">><br>
<b>Date: </b>Tuesday, 18 April 2023 at 11:08<br>
<b>To: </b>Chockalingam, Karthikeyan (STFC,DL,HC) <</span><a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank"><span style="font-size:12pt">karthikeyan.chockalingam@stfc.ac.uk</span></a><span style="font-size:12pt;color:black">><br>
<b>Cc: </b></span><a href="mailto:petsc-users@mcs.anl.gov" target="_blank"><span style="font-size:12pt">petsc-users@mcs.anl.gov</span></a><span style="font-size:12pt;color:black"> <</span><a href="mailto:petsc-users@mcs.anl.gov" target="_blank"><span style="font-size:12pt">petsc-users@mcs.anl.gov</span></a><span style="font-size:12pt;color:black">><br>
<b>Subject: </b>Re: [petsc-users] Setting up a matrix for Lagrange multiplier</span><u></u><u></u></p>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">On Tue, Apr 18, 2023 at 5:24 AM Karthikeyan Chockalingam - STFC UKRI via petsc-users <</span><a href="mailto:petsc-users@mcs.anl.gov" target="_blank"><span style="font-size:11pt">petsc-users@mcs.anl.gov</span></a><span style="font-size:11pt">>
wrote:</span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt;color:black">Hello,</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black">I'm solving a problem using the Lagrange multiplier, the matrix has the form</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black">K = [A P^T</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> P 0]</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt">I am familiar with constructing K using MATMPIAIJ. However, I would like to know if had [A], can I augment it with [P], [P^T] and [0] of type MATMPIAIJ?
Likewise for vectors as well.</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt">Can you please point me to the right resource, if it is a common operation in PETSc?</span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You can do this at least 2 ways:</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> 1) Assemble you submatrices directly into the larger matrix by constructing local-to-global maps for the emplacement. so that you do</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> not change your assembly code, except to change MatSetValues() to MatSetValuesLocal(). This is usually preferable.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> 2) Use MATNEST and VecNEST to put pointers to submatrices and subvectors directly in.</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Thanks,</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Matt</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><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:5pt 0cm 5pt 4.8pt">
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">Many thanks.</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt">Kind regards,</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt">Karthik.</span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><span style="font-size:11pt"><br clear="all">
</span><u></u><u></u></p>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<p class="MsoNormal"><span class="m_-968874379040273580m-4128250693350498318gmailsignatureprefix"><span style="font-size:11pt">--
</span></span><u></u><u></u></p>
<div>
<div>
<div>
<div>
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">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</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank"><span style="font-size:11pt">https://www.cse.buffalo.edu/~knepley/</span></a><u></u><u></u></p>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><span style="font-size:11pt"><br clear="all">
</span><u></u><u></u></p>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<p class="MsoNormal"><span class="m_-968874379040273580gmailsignatureprefix"><span style="font-size:11pt">--
</span></span><u></u><u></u></p>
<div>
<div>
<div>
<div>
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">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</span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> </span><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank"><span style="font-size:11pt">https://www.cse.buffalo.edu/~knepley/</span></a><u></u><u></u></p>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</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>