[petsc-users] question about MatSetLocalToGlobalMapping

Zhang, Hong hzhang at mcs.anl.gov
Fri Apr 21 10:29:21 CDT 2023


Pierre,
________________________________

1) Is there any hope to get PDIPDM to use a MatNest?

KKT matrix is indefinite and ill-conditioned, which must be solved using a direct matrix factorization method.
But you are using PCBJACOBI in the paper you attached?
In any case, there are many such systems, e.g., a discretization of Stokes equations, that can be solved with something else than a direct factorization.

The KKT in this paper is very special, each block represents a time step with loose and weak couplings. We tested larger and slightly varying physical parameters, PCZBJACOBI fails convergence while CHOLESKY encounters out of memory.
For the current implementation, we use MUMPS Cholesky as default. To use MatNest, what direct solver to use, SCHUR_FACTOR? I do not know how to get it work.
On the one hand, MatNest can efficiently convert to AIJ or SBAIJ if you want to stick to PCLU or PCCHOLESKY.
On the other hand, it allows to easily switch to PCFIELDSPLIT which can be used to solve saddle-point problems.

MatNest would be a natural way for these type of applications. PETSc has MatConvert_Nest_AIJ(), which could enable users to assemble their matrices in MatNest, and apply mumps/superlu direct solvers. I'm not concerned about the overhead of matrix conversion, because matrix factorization dominates computation in general.
2) Is this fixed https://lists.mcs.anl.gov/pipermail/petsc-dev/2020-September/026398.html ?
I cannot get users to transition away from Ipopt because of these two missing features.

The existing pdipm is the result of a MS student intern project. None of us involved are experts on the optimization solvers. We made a straightforward parallelization of Ipopt. It indeed needs further work, e.g., more features, better matrix storage, convergence criteria... To our knowledge, parallel pdipm is not available other than our pdipm.
Ipopt can use MUMPS and PARDISO internally, so it’s in some sense parallel (using shared memory).
Also, this is not a very potent selling point.
My users that are satisfied with Ipopt as a "non-parallel" black box don’t want to have to touch part of their code just to stick it in a parallel black box which is limited to the same kind of linear solver and which has severe limitations with respect to Hessian/Jacobian/constraint distributions.

We saw exiting 'parallel' pdipm papers, which conduct sequential computation and send KKT matrices to multiple processors, then they show scaling results on matrix computation only. Again, I wish someone can help to improve tao/pdipm. This is a useful solver.
We should improve our pdipm.
Hong

On 20 Apr 2023, at 5:47 PM, Zhang, Hong via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:

Karthik,
We built a KKT matrix in  TaoSNESJacobian_PDIPM() (see petsc/src/tao/constrained/impls/ipm/pdipm.c) which assembles several small matrices into a large KKT matrix in mpiaij format. You could take the same approach to insert P and P^T into your K.
FYI, I attached our paper.
Hong
________________________________

From: petsc-users <petsc-users-bounces at mcs.anl.gov<mailto:petsc-users-bounces at mcs.anl.gov>> on behalf of Matthew Knepley <knepley at gmail.com<mailto:knepley at gmail.com>>
Sent: Thursday, April 20, 2023 5:37 AM
To: Karthikeyan Chockalingam - STFC UKRI <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] question about MatSetLocalToGlobalMapping

On Thu, Apr 20, 2023 at 6:13 AM Karthikeyan Chockalingam - STFC UKRI via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:
Hello,

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)

K = [A P^T
       P   0]

Where K is of type MatMPIAIJ. I first constructed the top left [A] using MatSetValues().

Now, I would like to construct the bottom left [p] and top right [p^T] using MatSetValuesLocal().

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.

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

  (phi, phi)        --> A
  (phi, lambda) --> P^T
  (lambda, phi) --> P

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
maps. I do not like the row and col names since in P they reverse.

Q1) At what point should I declare MatSetLocalToGlobalMapping – is it just before I use MatSetValuesLocal()?

Okay, it is good you are asking this because my thinking was somewhat confused. I think the precise steps are:

  1) Create the large saddle point matrix K

    1a) We must call https://petsc.org/main/manualpages/Mat/MatSetLocalToGlobalMapping/ on it. In the simplest case, this just maps
           the local rows numbers [0, Nrows) to the global rows numbers [rowStart, rowStart + Nrows).

  2) To form each piece:

    2a) Extract that block using https://petsc.org/main/manualpages/Mat/MatGetLocalSubMatrix/

           This gives back a Mat object that you subsequently restore using https://petsc.org/main/manualpages/Mat/MatRestoreLocalSubMatrix/

     2b) Insert values using https://petsc.org/main/manualpages/Mat/MatSetValuesLocal/

            The local indices used for insertion here are indices relative to the block itself, and the L2G map for this matrix
            has been rewritten to insert into that block in the larger matrix. Thus this looks like just calling MatSetValuesLocal()
            on the smaller matrix block, but inserts correctly into the larger matrix.

Therefore, the code you write code in 2) could work equally well making the large matrix from 1), or independent smaller matrix blocks.

Does this make sense?

  Thanks,

     Matt

I will use MatSetLocalToGlobalMapping(K, row_mapping, column_mapping) to build the bottom left [P].


Q2) Can now I reset the mapping as MatSetLocalToGlobalMapping(K, column_mapping, row_mapping) to build the top right [P^T]?


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/>
<IET Generation Trans   Dist - 2022 - Sundermann - Parallel primal‐dual interior point method for the solution of dynamic.pdf>

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230421/6b34d158/attachment-0001.html>


More information about the petsc-users mailing list