<div dir="ltr"><div dir="ltr">On Mon, Apr 24, 2023 at 1:37 PM 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="msg7795183717874518080">
<div lang="EN-GB" style="overflow-wrap: break-word;">
<div class="m_7795183717874518080WordSection1">
<p class="MsoNormal"><span style="font-size:11pt">Great to know there is an example now.<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">Yes, -pc_fieldsplit_defect_saddle_point worked.<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">But I would like to understand what I did wrong (and why it didn’t work).</span></p></div></div></div></blockquote><div><br></div><div>Can you show the error? Or give something self-contained that I can run.</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="msg7795183717874518080"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_7795183717874518080WordSection1"><p class="MsoNormal"><span style="font-size:11pt">
<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">After reading many posted posts, I needed to create PCFieldSplitSetIS to split the fields.</span></p></div></div></div></blockquote><div><br></div><div>Yes. I would give NULL for the name here.</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="msg7795183717874518080"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_7795183717874518080WordSection1"><p class="MsoNormal"><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">I set the first N indices to the field phi and the next C indices (starting from N) to the field lambda.</span></p></div></div></div></blockquote><div><br></div><div>This looks fine. However, these are global indices, so this would not work in parallel. You would have to offset by</div><div>the first index on this process.</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="msg7795183717874518080"><div lang="EN-GB" style="overflow-wrap: break-word;"><div class="m_7795183717874518080WordSection1"><p class="MsoNormal"><span style="font-size:11pt">
<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">Please tell me what I did wrong and how I can fix the lines which are now commented?<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:rgb(57,0,160)"> PetscErrorCode ierr;<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(57,0,160)">KSP</span><span style="font-size:8.5pt;font-family:Menlo;color:black"> ksp;<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(57,0,160)">PC</span><span style="font-size:8.5pt;font-family:Menlo;color:black"> pc;<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"><u></u> <u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:rgb(108,54,169)"> KSPCreate(PETSC_COMM_WORLD, &ksp);<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(108,54,169)">KSPSetType</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(ksp,
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(100,56,32)">KSPGMRES</span><span style="font-size:8.5pt;font-family:Menlo;color:black">);<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(108,54,169)">KSPSetOperators</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(ksp,
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(50,109,116)">A</span><span style="font-size:8.5pt;font-family:Menlo;color:black">[level],
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(50,109,116)">A</span><span style="font-size:8.5pt;font-family:Menlo;color:black">[level]);<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"> ierr =
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(108,54,169)">KSPGetPC</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(ksp,&pc);</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(100,56,32)">CHKERRQ</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(ierr);<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"> ierr =
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(108,54,169)">PCSetType</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(pc,</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(100,56,32)">PCFIELDSPLIT</span><span style="font-size:8.5pt;font-family:Menlo;color:black">);</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(100,56,32)">CHKERRQ</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(ierr);<u></u><u></u></span></p>
<p class="m_7795183717874518080p1"><u></u> <u></u></p>
<p class="m_7795183717874518080p1"> ……..<u></u><u></u></p>
<p class="m_7795183717874518080p1"><u></u> <u></u></p>
<p class="MsoNormal" style="background:white"><b><span style="font-size:8.5pt;font-family:Menlo;color:rgb(155,35,147)"> for</span></b><span style="font-size:8.5pt;font-family:Menlo;color:black"> (</span><b><span style="font-size:8.5pt;font-family:Menlo;color:rgb(155,35,147)">int</span></b><span style="font-size:8.5pt;font-family:Menlo;color:black">
i = </span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(28,0,207)">0</span><span style="font-size:8.5pt;font-family:Menlo;color:black">; i < N; i++)<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"> vec_zero_field_index[i] = i;<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"> <u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(57,0,160)">PetscInt</span><span style="font-size:8.5pt;font-family:Menlo;color:black"> * vec_first_field_index =
</span><b><span style="font-size:8.5pt;font-family:Menlo;color:rgb(155,35,147)">NULL</span></b><span style="font-size:8.5pt;font-family:Menlo;color:black">;<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(100,56,32)">PetscMalloc</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(n *
</span><b><span style="font-size:8.5pt;font-family:Menlo;color:rgb(155,35,147)">sizeof</span></b><span style="font-size:8.5pt;font-family:Menlo;color:black">(</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(57,0,160)">PetscInt</span><span style="font-size:8.5pt;font-family:Menlo;color:black">),
&vec_first_field_index);<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"> <u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black">
</span><b><span style="font-size:8.5pt;font-family:Menlo;color:rgb(155,35,147)">for</span></b><span style="font-size:8.5pt;font-family:Menlo;color:black"> (</span><b><span style="font-size:8.5pt;font-family:Menlo;color:rgb(155,35,147)">int</span></b><span style="font-size:8.5pt;font-family:Menlo;color:black">
i = </span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(28,0,207)">0</span><span style="font-size:8.5pt;font-family:Menlo;color:black">; i < C; i++)<u></u><u></u></span></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:black"> vec_first_field_index[i] = N + i;<u></u><u></u></span></p>
<p class="m_7795183717874518080p1"><u></u> <u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080s1"> IS</span> zero_field_isx, first_field_isx;<u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span><span class="m_7795183717874518080s2">CHKERRQ</span>(<span class="m_7795183717874518080s3">ISCreateGeneral</span>(<span class="m_7795183717874518080s3">PETSC_COMM_WORLD</span>, N, vec_zero_field_index,
<span class="m_7795183717874518080s3">PETSC_COPY_VALUES</span>, &zero_field_isx));<u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span><span class="m_7795183717874518080s2">CHKERRQ</span>(<span class="m_7795183717874518080s3">ISCreateGeneral</span>(<span class="m_7795183717874518080s3">PETSC_COMM_WORLD</span>, C, vec_first_field_index,
<span class="m_7795183717874518080s3">PETSC_COPY_VALUES</span>, &first_field_isx));<u></u><u></u></p>
<p class="m_7795183717874518080p2"><span class="m_7795183717874518080apple-converted-space"> </span><u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = <span class="m_7795183717874518080s3">
PCFieldSplitSetIS</span>(pc,<span class="m_7795183717874518080s4">"0"</span>,zero_field_isx);<span class="m_7795183717874518080apple-converted-space">
</span><span class="m_7795183717874518080s2">CHKERRQ</span>(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = <span class="m_7795183717874518080s3">
PCFieldSplitSetIS</span>(pc,<span class="m_7795183717874518080s4">"1"</span>,first_field_isx);<span class="m_7795183717874518080apple-converted-space">
</span><span class="m_7795183717874518080s2">CHKERRQ</span>(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p2"><span class="m_7795183717874518080apple-converted-space"> </span><u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = <span class="m_7795183717874518080s3">
PetscOptionsSetValue</span>(<span class="m_7795183717874518080s5"><b>NULL</b></span>,<span class="m_7795183717874518080s4">"-ksp_type"</span>,
<span class="m_7795183717874518080s4">"fgmres"</span>); <span class="m_7795183717874518080s2">CHKERRQ</span>(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = <span class="m_7795183717874518080s3">
PetscOptionsSetValue</span>(<span class="m_7795183717874518080s5"><b>NULL</b></span>,<span class="m_7795183717874518080s4">"-pc_type"</span>,
<span class="m_7795183717874518080s4">"fieldsplit"</span>); <span class="m_7795183717874518080s2">CHKERRQ</span>(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p1"><u></u> <u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span><span class="m_7795183717874518080s6">/*</span><u></u><u></u></p>
<p class="m_7795183717874518080p3"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = PetscOptionsSetValue(NULL,"-fieldsplit_0_ksp_type", "gmres"); CHKERRQ(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p3"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = PetscOptionsSetValue(NULL,"-fieldsplit_0_pc_type", "jacobi"); CHKERRQ(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p3"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = PetscOptionsSetValue(NULL,"-fieldsplit_1_ksp_type", "preonly"); CHKERRQ(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p3"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = PetscOptionsSetValue(NULL,"-fieldsplit_1_pc_type", "lu"); CHKERRQ(ierr);*/<u></u><u></u></p>
<p class="m_7795183717874518080p3"><u></u> <u></u></p>
<p class="m_7795183717874518080p4"><span class="m_7795183717874518080apple-converted-space"> </span><span class="m_7795183717874518080s7">ierr =
</span><span class="m_7795183717874518080s3">PetscOptionsSetValue</span><span class="m_7795183717874518080s7">(</span><span class="m_7795183717874518080s5"><b>NULL</b></span><span class="m_7795183717874518080s7">,
</span>"-pc_fieldsplit_detect_saddle_point"<span class="m_7795183717874518080s7">, </span><span class="m_7795183717874518080s5"><b>NULL</b></span><span class="m_7795183717874518080s7">);</span><span class="m_7795183717874518080s2">CHKERRQ</span><span class="m_7795183717874518080s7">(ierr);</span><u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = <span class="m_7795183717874518080s3">
PetscOptionsSetValue</span>(<span class="m_7795183717874518080s5"><b>NULL</b></span>, <span class="m_7795183717874518080s4">
"-options_left"</span>, <span class="m_7795183717874518080s5"><b>NULL</b></span>);<span class="m_7795183717874518080s2">CHKERRQ</span>(ierr);<u></u><u></u></p>
<p class="m_7795183717874518080p1"><span class="m_7795183717874518080apple-converted-space"> </span>ierr = <span class="m_7795183717874518080s3">
KSPSetFromOptions</span>(ksp); <span class="m_7795183717874518080s2">CHKERRQ</span>(ierr);<u></u><u></u></p>
<p class="MsoNormal" style="background:white"><span style="font-size:8.5pt;font-family:Menlo;color:rgb(108,54,169)"> KSPSolve</span><span style="font-size:8.5pt;font-family:Menlo;color:black">(ksp,
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(50,109,116)">b</span><span style="font-size:8.5pt;font-family:Menlo;color:black">,
</span><span style="font-size:8.5pt;font-family:Menlo;color:rgb(50,109,116)">x</span><span style="font-size:8.5pt;font-family:Menlo;color:black">);<u></u><u></u></span></p>
<p class="m_7795183717874518080p1"><u></u> <u></u></p>
<p class="MsoNormal"><span style="font-size:11pt">Best,<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">Karthik.<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">Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
<b>Date: </b>Monday, 24 April 2023 at 17:57<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<u></u><u></u></span></p>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">On Mon, Apr 24, 2023 at 10:22 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></span></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>
<div>
<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><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> P 0]</span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt">and solved them using a direct solver.<u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">I modified your example to create either AIJ or Nest matrices, but use the same assembly code:<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <a href="https://gitlab.com/petsc/petsc/-/merge_requests/6368" target="_blank">https://gitlab.com/petsc/petsc/-/merge_requests/6368</a><u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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>
<div>
<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.<u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">With this particular matrix, you can use<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> -pc_fieldsplit_detect_saddle_point<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">and it will split it automatically.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Thanks,<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Matt<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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>
<div>
<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</span><span style="font-size:11pt"><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><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black"> </span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black"> ISLocalToGlobalMapping mapping;</span><span style="font-size:11pt"><u></u><u></u></span></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><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black"> </span><span style="font-size:11pt"><u></u><u></u></span></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><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black"> </span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black">Thank you.</span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black"> </span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black">Kind regards,</span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal" style="background:white">
<span style="font-size:11pt;color:black">Karthik.</span><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>
<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">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><span style="font-size:11pt"><u></u><u></u></span></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 <<a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank">karthikeyan.chockalingam@stfc.ac.uk</a>>
wrote:<u></u><u></u></span></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_7795183717874518080m-968874379040273580m-4128250693350498318p1">Thank you for your response. I spend some time understanding how<u></u><u></u></p>
<p class="m_7795183717874518080m-968874379040273580m-4128250693350498318p1">MatSetValuesLocal and ISLocalToGlobalMappingCreate work.<u></u><u></u></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You can look at SNES ex28 where we do this with DMCOMPOSITE.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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?<u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black">K = [A P^T</span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> P 0]</span><span style="font-size:11pt"><u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">I assume MPIAIJ since IS is only used for Neumann-Neumann decompositions.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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.<u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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?<u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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<u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You wouldn't in this case.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Thanks,<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Matt<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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,<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>
<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">Matthew Knepley <</span><span style="font-size:11pt"><a href="mailto:knepley@gmail.com" target="_blank"><span style="font-size:12pt">knepley@gmail.com</span></a></span><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><span style="font-size:11pt"><a href="mailto:karthikeyan.chockalingam@stfc.ac.uk" target="_blank"><span style="font-size:12pt">karthikeyan.chockalingam@stfc.ac.uk</span></a></span><span style="font-size:12pt;color:black">><br>
<b>Cc: </b></span><span style="font-size:11pt"><a href="mailto:petsc-users@mcs.anl.gov" target="_blank"><span style="font-size:12pt">petsc-users@mcs.anl.gov</span></a></span><span style="font-size:12pt;color:black"> <</span><span style="font-size:11pt"><a href="mailto:petsc-users@mcs.anl.gov" target="_blank"><span style="font-size:12pt">petsc-users@mcs.anl.gov</span></a></span><span style="font-size:12pt;color:black">><br>
<b>Subject: </b>Re: [petsc-users] Setting up a matrix for Lagrange multiplier</span><span style="font-size:11pt"><u></u><u></u></span></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 <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>>
wrote:<u></u><u></u></span></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><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> </span><span style="font-size:11pt"><u></u><u></u></span></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><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> </span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black">K = [A P^T</span><span style="font-size:11pt"><u></u><u></u></span></p>
<p class="MsoNormal"><span style="font-size:11pt;color:black"> P 0]</span><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>
<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.<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">Can you please point me to the right resource, if it is a common operation in PETSc?<u></u><u></u></span></p>
</div>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt">You can do this at least 2 ways:<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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<u></u><u></u></span></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.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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.<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Thanks,<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> Matt<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></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.<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">Kind 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>
<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>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><span style="font-size:11pt"><br clear="all">
<u></u><u></u></span></p>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<p class="MsoNormal"><span class="m_7795183717874518080m-968874379040273580m-4128250693350498318gmailsignatureprefix"><span style="font-size:11pt">--
</span></span><span style="font-size:11pt"><u></u><u></u></span></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<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><u></u><u></u></span></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">
<u></u><u></u></span></p>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<p class="MsoNormal"><span class="m_7795183717874518080m-968874379040273580gmailsignatureprefix"><span style="font-size:11pt">--
</span></span><span style="font-size:11pt"><u></u><u></u></span></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<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"> <u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><u></u><u></u></span></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">
<u></u><u></u></span></p>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<p class="MsoNormal"><span class="m_7795183717874518080gmailsignatureprefix"><span style="font-size:11pt">--
</span></span><span style="font-size:11pt"><u></u><u></u></span></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<u></u><u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><u></u> <u></u></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:11pt"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><u></u><u></u></span></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>