[petsc-users] Setting up a matrix for Lagrange multiplier

Karthikeyan Chockalingam - STFC UKRI karthikeyan.chockalingam at stfc.ac.uk
Mon Apr 24 09:21:57 CDT 2023


Hello,

I was able to construct the below K matrix (using submatrices P and P^T), which is of type MATAIJ
K = [A P^T
       P   0]
and solved them using a direct solver.
However, I was reading online that this is a saddle point problem and I should be employing PCFIELDSPLIT.
Since I have one monolithic matrix K, I was not sure how to split the fields.

Best regards,
Karthik.


From: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalingam at stfc.ac.uk>
Date: Wednesday, 19 April 2023 at 17:52
To: Matthew Knepley <knepley at gmail.com>
Cc: petsc-users at mcs.anl.gov <petsc-users at mcs.anl.gov>
Subject: Re: [petsc-users] Setting up a matrix for Lagrange multiplier
 I have declared the mapping

  ISLocalToGlobalMapping mapping;
    ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, n, nindices, PETSC_COPY_VALUES, &mapping);

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?

Thank you.

Kind regards,
Karthik.


From: Matthew Knepley <knepley at gmail.com>
Date: Tuesday, 18 April 2023 at 16:21
To: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalingam at stfc.ac.uk>
Cc: petsc-users at mcs.anl.gov <petsc-users at mcs.anl.gov>
Subject: Re: [petsc-users] Setting up a matrix for Lagrange multiplier
On Tue, Apr 18, 2023 at 11:16 AM Karthikeyan Chockalingam - STFC UKRI <karthikeyan.chockalingam at stfc.ac.uk<mailto:karthikeyan.chockalingam at stfc.ac.uk>> wrote:

Thank you for your response. I spend some time understanding how

MatSetValuesLocal and ISLocalToGlobalMappingCreate work.

You can look at SNES ex28 where we do this with DMCOMPOSITE.

Q1) Will the matrix K be of type MATMPIAIJ or MATIS?
K = [A P^T
       P   0]

I assume MPIAIJ since IS is only used for Neumann-Neumann decompositions.

Q2) Can I use both MatSetValues() to MatSetValuesLocal() to populate K? Since I have already used MatSetValues() to construct A.

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.

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?

You have a monolithic K, so that you can use sparse direct solvers to check things. THis is impossible with separate storage.

Q4) I probably don’t have to construct an independent P matrix

You wouldn't in this case.

  Thanks,

     Matt

Best regards,
Karthik.



From: Matthew Knepley <knepley at gmail.com<mailto:knepley at gmail.com>>
Date: Tuesday, 18 April 2023 at 11:08
To: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalingam at stfc.ac.uk<mailto:karthikeyan.chockalingam at stfc.ac.uk>>
Cc: petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov> <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>>
Subject: Re: [petsc-users] Setting up a matrix for Lagrange multiplier
On Tue, Apr 18, 2023 at 5:24 AM Karthikeyan Chockalingam - STFC UKRI via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:
Hello,

I'm solving a problem using the Lagrange multiplier, the matrix has the form

K = [A P^T
       P   0]

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.

Can you please point me to the right resource, if it is a common operation in PETSc?

You can do this at least 2 ways:

  1) Assemble you submatrices directly into the larger matrix by constructing local-to-global maps for the emplacement. so that you do
      not change your assembly code, except to change MatSetValues() to MatSetValuesLocal(). This is usually preferable.

  2) Use MATNEST and VecNEST to put pointers to submatrices and subvectors directly in.

  Thanks,

     Matt

Many thanks.

Kind regards,
Karthik.






--
What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/<http://www.cse.buffalo.edu/~knepley/>


--
What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/<http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230424/bfcafe81/attachment-0001.html>


More information about the petsc-users mailing list