xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
567602552SJunchao Zhang #include <petscpkg_version.h>
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
87ee00b23SStefano Zampini #include <../src/mat/impls/sell/mpi/mpisell.h>
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
122907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
132907cef9SHong Zhang     #include <cmumps_c.h>
142907cef9SHong Zhang   #else
15c6db04a5SJed Brown     #include <zmumps_c.h>
162907cef9SHong Zhang   #endif
172907cef9SHong Zhang #else
182907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
192907cef9SHong Zhang     #include <smumps_c.h>
20397b6df1SKris Buschelman   #else
21c6db04a5SJed Brown     #include <dmumps_c.h>
22397b6df1SKris Buschelman   #endif
232907cef9SHong Zhang #endif
24397b6df1SKris Buschelman EXTERN_C_END
25397b6df1SKris Buschelman #define JOB_INIT         -1
26413bcc21SPierre Jolivet #define JOB_NULL         0
273d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
283d472b54SHong Zhang #define JOB_FACTNUMERIC  2
293d472b54SHong Zhang #define JOB_SOLVE        3
30397b6df1SKris Buschelman #define JOB_END          -2
313d472b54SHong Zhang 
322907cef9SHong Zhang /* calls to MUMPS */
332907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
342907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
353ab56b82SJunchao Zhang     #define MUMPS_c cmumps_c
362907cef9SHong Zhang   #else
373ab56b82SJunchao Zhang     #define MUMPS_c zmumps_c
382907cef9SHong Zhang   #endif
392907cef9SHong Zhang #else
402907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
413ab56b82SJunchao Zhang     #define MUMPS_c smumps_c
422907cef9SHong Zhang   #else
433ab56b82SJunchao Zhang     #define MUMPS_c dmumps_c
442907cef9SHong Zhang   #endif
452907cef9SHong Zhang #endif
462907cef9SHong Zhang 
47a6053eceSJunchao Zhang /* MUMPS uses MUMPS_INT for nonzero indices such as irn/jcn, irn_loc/jcn_loc and uses int64_t for
48a6053eceSJunchao Zhang    number of nonzeros such as nnz, nnz_loc. We typedef MUMPS_INT to PetscMUMPSInt to follow the
49a6053eceSJunchao Zhang    naming convention in PetscMPIInt, PetscBLASInt etc.
50a6053eceSJunchao Zhang */
51a6053eceSJunchao Zhang typedef MUMPS_INT PetscMUMPSInt;
52a6053eceSJunchao Zhang 
5367602552SJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_GE(5, 3, 0)
5467602552SJunchao Zhang   #if defined(MUMPS_INTSIZE64) /* MUMPS_INTSIZE64 is in MUMPS headers if it is built in full 64-bit mode, therefore the macro is more reliable */
55a6053eceSJunchao Zhang     #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
5667602552SJunchao Zhang   #endif
57a6053eceSJunchao Zhang #else
5867602552SJunchao Zhang   #if defined(INTSIZE64) /* INTSIZE64 is a command line macro one used to build MUMPS in full 64-bit mode */
5967602552SJunchao Zhang     #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
6067602552SJunchao Zhang   #endif
6167602552SJunchao Zhang #endif
6267602552SJunchao Zhang 
63a6053eceSJunchao Zhang #define MPIU_MUMPSINT       MPI_INT
64a6053eceSJunchao Zhang #define PETSC_MUMPS_INT_MAX 2147483647
65a6053eceSJunchao Zhang #define PETSC_MUMPS_INT_MIN -2147483648
66a6053eceSJunchao Zhang 
67a6053eceSJunchao Zhang /* Cast PetscInt to PetscMUMPSInt. Usually there is no overflow since <a> is row/col indices or some small integers*/
68d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscMUMPSIntCast(PetscInt a, PetscMUMPSInt *b)
69d71ae5a4SJacob Faibussowitsch {
70a6053eceSJunchao Zhang   PetscFunctionBegin;
71ece88022SPierre Jolivet #if PetscDefined(USE_64BIT_INDICES)
722c71b3e2SJacob Faibussowitsch   PetscAssert(a <= PETSC_MUMPS_INT_MAX && a >= PETSC_MUMPS_INT_MIN, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "PetscInt too long for PetscMUMPSInt");
73ece88022SPierre Jolivet #endif
74a6053eceSJunchao Zhang   *b = (PetscMUMPSInt)(a);
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76a6053eceSJunchao Zhang }
77a6053eceSJunchao Zhang 
78a6053eceSJunchao Zhang /* Put these utility routines here since they are only used in this file */
79d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscOptionsMUMPSInt_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscMUMPSInt currentvalue, PetscMUMPSInt *value, PetscBool *set, PetscMUMPSInt lb, PetscMUMPSInt ub)
80d71ae5a4SJacob Faibussowitsch {
81a6053eceSJunchao Zhang   PetscInt  myval;
82a6053eceSJunchao Zhang   PetscBool myset;
83a6053eceSJunchao Zhang   PetscFunctionBegin;
84a6053eceSJunchao Zhang   /* PetscInt's size should be always >= PetscMUMPSInt's. It is safe to call PetscOptionsInt_Private to read a PetscMUMPSInt */
859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, (PetscInt)currentvalue, &myval, &myset, lb, ub));
869566063dSJacob Faibussowitsch   if (myset) PetscCall(PetscMUMPSIntCast(myval, value));
87a6053eceSJunchao Zhang   if (set) *set = myset;
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89a6053eceSJunchao Zhang }
90a6053eceSJunchao Zhang #define PetscOptionsMUMPSInt(a, b, c, d, e, f) PetscOptionsMUMPSInt_Private(PetscOptionsObject, a, b, c, d, e, f, PETSC_MUMPS_INT_MIN, PETSC_MUMPS_INT_MAX)
91a6053eceSJunchao Zhang 
92217d3b1eSJunchao Zhang /* if using PETSc OpenMP support, we only call MUMPS on master ranks. Before/after the call, we change/restore CPUs the master ranks can run on */
933ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
943ab56b82SJunchao Zhang   #define PetscMUMPS_c(mumps) \
953ab56b82SJunchao Zhang     do { \
963ab56b82SJunchao Zhang       if (mumps->use_petsc_omp_support) { \
973ab56b82SJunchao Zhang         if (mumps->is_omp_master) { \
989566063dSJacob Faibussowitsch           PetscCall(PetscOmpCtrlOmpRegionOnMasterBegin(mumps->omp_ctrl)); \
9914ffdc6fSStefano Zampini           PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
10014ffdc6fSStefano Zampini           PetscStackCallExternalVoid(PetscStringize(MUMPS_c), MUMPS_c(&mumps->id)); \
10114ffdc6fSStefano Zampini           PetscCall(PetscFPTrapPop()); \
1029566063dSJacob Faibussowitsch           PetscCall(PetscOmpCtrlOmpRegionOnMasterEnd(mumps->omp_ctrl)); \
1033ab56b82SJunchao Zhang         } \
1049566063dSJacob Faibussowitsch         PetscCall(PetscOmpCtrlBarrier(mumps->omp_ctrl)); \
105c3714a1dSJunchao Zhang         /* Global info is same on all processes so we Bcast it within omp_comm. Local info is specific      \
106c3714a1dSJunchao Zhang          to processes, so we only Bcast info[1], an error code and leave others (since they do not have   \
107c3714a1dSJunchao Zhang          an easy translation between omp_comm and petsc_comm). See MUMPS-5.1.2 manual p82.                   \
108c3714a1dSJunchao Zhang          omp_comm is a small shared memory communicator, hence doing multiple Bcast as shown below is OK. \
109c3714a1dSJunchao Zhang       */ \
1109566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Bcast(mumps->id.infog, 40, MPIU_MUMPSINT, 0, mumps->omp_comm)); \
1119566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Bcast(mumps->id.rinfog, 20, MPIU_REAL, 0, mumps->omp_comm)); \
1129566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Bcast(mumps->id.info, 1, MPIU_MUMPSINT, 0, mumps->omp_comm)); \
1133ab56b82SJunchao Zhang       } else { \
11414ffdc6fSStefano Zampini         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
11514ffdc6fSStefano Zampini         PetscStackCallExternalVoid(PetscStringize(MUMPS_c), MUMPS_c(&mumps->id)); \
11614ffdc6fSStefano Zampini         PetscCall(PetscFPTrapPop()); \
1173ab56b82SJunchao Zhang       } \
1183ab56b82SJunchao Zhang     } while (0)
1193ab56b82SJunchao Zhang #else
1203ab56b82SJunchao Zhang   #define PetscMUMPS_c(mumps) \
121d71ae5a4SJacob Faibussowitsch     do { \
12214ffdc6fSStefano Zampini       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
12314ffdc6fSStefano Zampini       PetscStackCallExternalVoid(PetscStringize(MUMPS_c), MUMPS_c(&mumps->id)); \
12414ffdc6fSStefano Zampini       PetscCall(PetscFPTrapPop()); \
125d71ae5a4SJacob Faibussowitsch     } while (0)
1263ab56b82SJunchao Zhang #endif
1273ab56b82SJunchao Zhang 
128940cd9d6SSatish Balay /* declare MumpsScalar */
129940cd9d6SSatish Balay #if defined(PETSC_USE_COMPLEX)
130940cd9d6SSatish Balay   #if defined(PETSC_USE_REAL_SINGLE)
131940cd9d6SSatish Balay     #define MumpsScalar mumps_complex
132940cd9d6SSatish Balay   #else
133940cd9d6SSatish Balay     #define MumpsScalar mumps_double_complex
134940cd9d6SSatish Balay   #endif
135940cd9d6SSatish Balay #else
136940cd9d6SSatish Balay   #define MumpsScalar PetscScalar
137940cd9d6SSatish Balay #endif
1383d472b54SHong Zhang 
139397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
140397b6df1SKris Buschelman #define ICNTL(I)  icntl[(I)-1]
141397b6df1SKris Buschelman #define CNTL(I)   cntl[(I)-1]
142397b6df1SKris Buschelman #define INFOG(I)  infog[(I)-1]
143a7aca84bSHong Zhang #define INFO(I)   info[(I)-1]
144397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
145adc1d99fSHong Zhang #define RINFO(I)  rinfo[(I)-1]
146397b6df1SKris Buschelman 
147a6053eceSJunchao Zhang typedef struct Mat_MUMPS Mat_MUMPS;
148a6053eceSJunchao Zhang struct Mat_MUMPS {
149397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
1502907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
1512907cef9SHong Zhang   CMUMPS_STRUC_C id;
1522907cef9SHong Zhang   #else
153397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
1542907cef9SHong Zhang   #endif
1552907cef9SHong Zhang #else
1562907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
1572907cef9SHong Zhang   SMUMPS_STRUC_C id;
158397b6df1SKris Buschelman   #else
159397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
160397b6df1SKris Buschelman   #endif
1612907cef9SHong Zhang #endif
1622907cef9SHong Zhang 
163397b6df1SKris Buschelman   MatStructure   matstruc;
1642d4298aeSJunchao Zhang   PetscMPIInt    myid, petsc_size;
165a6053eceSJunchao Zhang   PetscMUMPSInt *irn, *jcn;       /* the (i,j,v) triplets passed to mumps. */
166a6053eceSJunchao Zhang   PetscScalar   *val, *val_alloc; /* For some matrices, we can directly access their data array without a buffer. For others, we need a buffer. So comes val_alloc. */
167a6053eceSJunchao Zhang   PetscInt64     nnz;             /* number of nonzeros. The type is called selective 64-bit in mumps */
168a6053eceSJunchao Zhang   PetscMUMPSInt  sym;
1692d4298aeSJunchao Zhang   MPI_Comm       mumps_comm;
170413bcc21SPierre Jolivet   PetscMUMPSInt *ICNTL_pre;
171413bcc21SPierre Jolivet   PetscReal     *CNTL_pre;
172a6053eceSJunchao Zhang   PetscMUMPSInt  ICNTL9_pre;         /* check if ICNTL(9) is changed from previous MatSolve */
173801fbe65SHong Zhang   VecScatter     scat_rhs, scat_sol; /* used by MatSolve() */
17425aac85cSJunchao Zhang   PetscMUMPSInt  ICNTL20;            /* use centralized (0) or distributed (10) dense RHS */
17567602552SJunchao Zhang   PetscMUMPSInt  lrhs_loc, nloc_rhs, *irhs_loc;
17667602552SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
17767602552SJunchao Zhang   PetscInt    *rhs_nrow, max_nrhs;
17867602552SJunchao Zhang   PetscMPIInt *rhs_recvcounts, *rhs_disps;
17967602552SJunchao Zhang   PetscScalar *rhs_loc, *rhs_recvbuf;
18067602552SJunchao Zhang #endif
181801fbe65SHong Zhang   Vec            b_seq, x_seq;
182a6053eceSJunchao Zhang   PetscInt       ninfo, *info; /* which INFO to display */
183b5fa320bSStefano Zampini   PetscInt       sizeredrhs;
18459ac8732SStefano Zampini   PetscScalar   *schur_sol;
18559ac8732SStefano Zampini   PetscInt       schur_sizesol;
186a6053eceSJunchao Zhang   PetscMUMPSInt *ia_alloc, *ja_alloc; /* work arrays used for the CSR struct for sparse rhs */
187a6053eceSJunchao Zhang   PetscInt64     cur_ilen, cur_jlen;  /* current len of ia_alloc[], ja_alloc[] */
188a6053eceSJunchao Zhang   PetscErrorCode (*ConvertToTriples)(Mat, PetscInt, MatReuse, Mat_MUMPS *);
1892205254eSKarl Rupp 
190a6053eceSJunchao Zhang   /* stuff used by petsc/mumps OpenMP support*/
1913ab56b82SJunchao Zhang   PetscBool    use_petsc_omp_support;
192da81f932SPierre Jolivet   PetscOmpCtrl omp_ctrl;             /* an OpenMP controller that blocked processes will release their CPU (MPI_Barrier does not have this guarantee) */
1933ab56b82SJunchao Zhang   MPI_Comm     petsc_comm, omp_comm; /* petsc_comm is petsc matrix's comm */
194a6053eceSJunchao Zhang   PetscInt64  *recvcount;            /* a collection of nnz on omp_master */
195a6053eceSJunchao Zhang   PetscMPIInt  tag, omp_comm_size;
1963ab56b82SJunchao Zhang   PetscBool    is_omp_master; /* is this rank the master of omp_comm */
197a6053eceSJunchao Zhang   MPI_Request *reqs;
198a6053eceSJunchao Zhang };
1993ab56b82SJunchao Zhang 
200a6053eceSJunchao Zhang /* Cast a 1-based CSR represented by (nrow, ia, ja) of type PetscInt to a CSR of type PetscMUMPSInt.
201a6053eceSJunchao Zhang    Here, nrow is number of rows, ia[] is row pointer and ja[] is column indices.
202a6053eceSJunchao Zhang  */
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscMUMPSIntCSRCast(Mat_MUMPS *mumps, PetscInt nrow, PetscInt *ia, PetscInt *ja, PetscMUMPSInt **ia_mumps, PetscMUMPSInt **ja_mumps, PetscMUMPSInt *nnz_mumps)
204d71ae5a4SJacob Faibussowitsch {
205a6053eceSJunchao Zhang   PetscInt nnz = ia[nrow] - 1; /* mumps uses 1-based indices. Uses PetscInt instead of PetscInt64 since mumps only uses PetscMUMPSInt for rhs */
206f0c56d0fSKris Buschelman 
207a6053eceSJunchao Zhang   PetscFunctionBegin;
208a6053eceSJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
209a6053eceSJunchao Zhang   {
210a6053eceSJunchao Zhang     PetscInt i;
211a6053eceSJunchao Zhang     if (nrow + 1 > mumps->cur_ilen) { /* realloc ia_alloc/ja_alloc to fit ia/ja */
2129566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->ia_alloc));
2139566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrow + 1, &mumps->ia_alloc));
214a6053eceSJunchao Zhang       mumps->cur_ilen = nrow + 1;
215a6053eceSJunchao Zhang     }
216a6053eceSJunchao Zhang     if (nnz > mumps->cur_jlen) {
2179566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->ja_alloc));
2189566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nnz, &mumps->ja_alloc));
219a6053eceSJunchao Zhang       mumps->cur_jlen = nnz;
220a6053eceSJunchao Zhang     }
2219566063dSJacob Faibussowitsch     for (i = 0; i < nrow + 1; i++) PetscCall(PetscMUMPSIntCast(ia[i], &(mumps->ia_alloc[i])));
2229566063dSJacob Faibussowitsch     for (i = 0; i < nnz; i++) PetscCall(PetscMUMPSIntCast(ja[i], &(mumps->ja_alloc[i])));
223a6053eceSJunchao Zhang     *ia_mumps = mumps->ia_alloc;
224a6053eceSJunchao Zhang     *ja_mumps = mumps->ja_alloc;
225a6053eceSJunchao Zhang   }
226a6053eceSJunchao Zhang #else
227a6053eceSJunchao Zhang   *ia_mumps          = ia;
228a6053eceSJunchao Zhang   *ja_mumps          = ja;
229a6053eceSJunchao Zhang #endif
2309566063dSJacob Faibussowitsch   PetscCall(PetscMUMPSIntCast(nnz, nnz_mumps));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232a6053eceSJunchao Zhang }
233b24902e0SBarry Smith 
234d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS *mumps)
235d71ae5a4SJacob Faibussowitsch {
236b5fa320bSStefano Zampini   PetscFunctionBegin;
2379566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.listvar_schur));
2389566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.redrhs));
2399566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->schur_sol));
24059ac8732SStefano Zampini   mumps->id.size_schur = 0;
241b3cb21ddSStefano Zampini   mumps->id.schur_lld  = 0;
24259ac8732SStefano Zampini   mumps->id.ICNTL(19)  = 0;
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24459ac8732SStefano Zampini }
24559ac8732SStefano Zampini 
246b3cb21ddSStefano Zampini /* solve with rhs in mumps->id.redrhs and return in the same location */
247d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
248d71ae5a4SJacob Faibussowitsch {
249b3cb21ddSStefano Zampini   Mat_MUMPS           *mumps = (Mat_MUMPS *)F->data;
250b3cb21ddSStefano Zampini   Mat                  S, B, X;
251b3cb21ddSStefano Zampini   MatFactorSchurStatus schurstatus;
252b3cb21ddSStefano Zampini   PetscInt             sizesol;
25359ac8732SStefano Zampini 
25459ac8732SStefano Zampini   PetscFunctionBegin;
2559566063dSJacob Faibussowitsch   PetscCall(MatFactorFactorizeSchurComplement(F));
2569566063dSJacob Faibussowitsch   PetscCall(MatFactorGetSchurComplement(F, &S, &schurstatus));
2579566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mumps->id.size_schur, mumps->id.nrhs, (PetscScalar *)mumps->id.redrhs, &B));
2589566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, ((PetscObject)S)->type_name));
259a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2609566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(B, S->boundtocpu));
261a3d589ffSStefano Zampini #endif
262b3cb21ddSStefano Zampini   switch (schurstatus) {
263d71ae5a4SJacob Faibussowitsch   case MAT_FACTOR_SCHUR_FACTORED:
264d71ae5a4SJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mumps->id.size_schur, mumps->id.nrhs, (PetscScalar *)mumps->id.redrhs, &X));
265d71ae5a4SJacob Faibussowitsch     PetscCall(MatSetType(X, ((PetscObject)S)->type_name));
266a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2679566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(X, S->boundtocpu));
268a3d589ffSStefano Zampini #endif
269b3cb21ddSStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
2709566063dSJacob Faibussowitsch       PetscCall(MatMatSolveTranspose(S, B, X));
27159ac8732SStefano Zampini     } else {
2729566063dSJacob Faibussowitsch       PetscCall(MatMatSolve(S, B, X));
27359ac8732SStefano Zampini     }
274b3cb21ddSStefano Zampini     break;
275b3cb21ddSStefano Zampini   case MAT_FACTOR_SCHUR_INVERTED:
276b3cb21ddSStefano Zampini     sizesol = mumps->id.nrhs * mumps->id.size_schur;
27759ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
2789566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->schur_sol));
2799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sizesol, &mumps->schur_sol));
28059ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
281b5fa320bSStefano Zampini     }
2829566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mumps->id.size_schur, mumps->id.nrhs, mumps->schur_sol, &X));
2839566063dSJacob Faibussowitsch     PetscCall(MatSetType(X, ((PetscObject)S)->type_name));
284a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2859566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(X, S->boundtocpu));
286a3d589ffSStefano Zampini #endif
2879566063dSJacob Faibussowitsch     PetscCall(MatProductCreateWithMat(S, B, NULL, X));
28859ac8732SStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
2899566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(X, MATPRODUCT_AtB));
290b5fa320bSStefano Zampini     } else {
2919566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(X, MATPRODUCT_AB));
292b5fa320bSStefano Zampini     }
2939566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(X));
2949566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(X));
2959566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(X));
2964417c5e8SHong Zhang 
2979566063dSJacob Faibussowitsch     PetscCall(MatCopy(X, B, SAME_NONZERO_PATTERN));
298b3cb21ddSStefano Zampini     break;
299d71ae5a4SJacob Faibussowitsch   default:
300d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status);
30159ac8732SStefano Zampini   }
3029566063dSJacob Faibussowitsch   PetscCall(MatFactorRestoreSchurComplement(F, &S, schurstatus));
3039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
3049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306b5fa320bSStefano Zampini }
307b5fa320bSStefano Zampini 
308d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
309d71ae5a4SJacob Faibussowitsch {
310b3cb21ddSStefano Zampini   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
311b5fa320bSStefano Zampini 
312b5fa320bSStefano Zampini   PetscFunctionBegin;
313b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
3143ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
315b5fa320bSStefano Zampini   }
316b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
317b5fa320bSStefano Zampini     PetscInt sizeredrhs = mumps->id.nrhs * mumps->id.size_schur;
318b5fa320bSStefano Zampini     /* allocate MUMPS internal array to store reduced right-hand sides */
319b5fa320bSStefano Zampini     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
3209566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->id.redrhs));
321b5fa320bSStefano Zampini       mumps->id.lredrhs = mumps->id.size_schur;
3229566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mumps->id.nrhs * mumps->id.lredrhs, &mumps->id.redrhs));
323b5fa320bSStefano Zampini       mumps->sizeredrhs = mumps->id.nrhs * mumps->id.lredrhs;
324b5fa320bSStefano Zampini     }
325b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 1; /* condensation phase */
326b5fa320bSStefano Zampini   } else {                   /* prepare for the expansion step */
327b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
3289566063dSJacob Faibussowitsch     PetscCall(MatMumpsSolveSchur_Private(F));
329b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
3303ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
33108401ef6SPierre Jolivet     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in solve phase: INFOG(1)=%d", mumps->id.INFOG(1));
332b5fa320bSStefano Zampini     /* restore defaults */
333b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
334d3d598ffSStefano Zampini     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
335d3d598ffSStefano Zampini     if (mumps->id.nrhs > 1) {
3369566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->id.redrhs));
337d3d598ffSStefano Zampini       mumps->id.lredrhs = 0;
338d3d598ffSStefano Zampini       mumps->sizeredrhs = 0;
339d3d598ffSStefano Zampini     }
340b5fa320bSStefano Zampini   }
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
342b5fa320bSStefano Zampini }
343b5fa320bSStefano Zampini 
344397b6df1SKris Buschelman /*
345d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
346d341cd04SHong Zhang 
347397b6df1SKris Buschelman   input:
34875480915SPierre Jolivet     A       - matrix in aij,baij or sbaij format
349397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
350bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
351bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
352397b6df1SKris Buschelman   output:
353397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
354397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
355eb9baa12SBarry Smith 
356eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
3577ee00b23SStefano Zampini   freed with PetscFree(mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
358eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
359eb9baa12SBarry Smith 
360397b6df1SKris Buschelman  */
36116ebf90aSShri Abhyankar 
362d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
363d71ae5a4SJacob Faibussowitsch {
364a3d589ffSStefano Zampini   const PetscScalar *av;
365185f6596SHong Zhang   const PetscInt    *ai, *aj, *ajj, M = A->rmap->n;
366a6053eceSJunchao Zhang   PetscInt64         nz, rnz, i, j, k;
367a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
36816ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa = (Mat_SeqAIJ *)A->data;
369397b6df1SKris Buschelman 
370397b6df1SKris Buschelman   PetscFunctionBegin;
3719566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(A, &av));
372a6053eceSJunchao Zhang   mumps->val = (PetscScalar *)av;
373bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
3742205254eSKarl Rupp     nz = aa->nz;
3752205254eSKarl Rupp     ai = aa->i;
3762205254eSKarl Rupp     aj = aa->j;
3779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
378a6053eceSJunchao Zhang     for (i = k = 0; i < M; i++) {
37916ebf90aSShri Abhyankar       rnz = ai[i + 1] - ai[i];
38067877ebaSShri Abhyankar       ajj = aj + ai[i];
38167877ebaSShri Abhyankar       for (j = 0; j < rnz; j++) {
3829566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(i + shift, &row[k]));
3839566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[k]));
384a6053eceSJunchao Zhang         k++;
38516ebf90aSShri Abhyankar       }
38616ebf90aSShri Abhyankar     }
387a6053eceSJunchao Zhang     mumps->irn = row;
388a6053eceSJunchao Zhang     mumps->jcn = col;
389a6053eceSJunchao Zhang     mumps->nnz = nz;
39016ebf90aSShri Abhyankar   }
3919566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39316ebf90aSShri Abhyankar }
394397b6df1SKris Buschelman 
395d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
396d71ae5a4SJacob Faibussowitsch {
397a6053eceSJunchao Zhang   PetscInt64     nz, i, j, k, r;
3987ee00b23SStefano Zampini   Mat_SeqSELL   *a = (Mat_SeqSELL *)A->data;
399a6053eceSJunchao Zhang   PetscMUMPSInt *row, *col;
4007ee00b23SStefano Zampini 
4017ee00b23SStefano Zampini   PetscFunctionBegin;
402a6053eceSJunchao Zhang   mumps->val = a->val;
4037ee00b23SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
4047ee00b23SStefano Zampini     nz = a->sliidx[a->totalslices];
4059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
406a6053eceSJunchao Zhang     for (i = k = 0; i < a->totalslices; i++) {
40748a46eb9SPierre Jolivet       for (j = a->sliidx[i], r = 0; j < a->sliidx[i + 1]; j++, r = ((r + 1) & 0x07)) PetscCall(PetscMUMPSIntCast(8 * i + r + shift, &row[k++]));
4087ee00b23SStefano Zampini     }
4099566063dSJacob Faibussowitsch     for (i = 0; i < nz; i++) PetscCall(PetscMUMPSIntCast(a->colidx[i] + shift, &col[i]));
410a6053eceSJunchao Zhang     mumps->irn = row;
411a6053eceSJunchao Zhang     mumps->jcn = col;
412a6053eceSJunchao Zhang     mumps->nnz = nz;
4137ee00b23SStefano Zampini   }
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4157ee00b23SStefano Zampini }
4167ee00b23SStefano Zampini 
417d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
418d71ae5a4SJacob Faibussowitsch {
41967877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa = (Mat_SeqBAIJ *)A->data;
42033d57670SJed Brown   const PetscInt *ai, *aj, *ajj, bs2 = aa->bs2;
421a6053eceSJunchao Zhang   PetscInt64      M, nz, idx = 0, rnz, i, j, k, m;
422a6053eceSJunchao Zhang   PetscInt        bs;
423a6053eceSJunchao Zhang   PetscMUMPSInt  *row, *col;
42467877ebaSShri Abhyankar 
42567877ebaSShri Abhyankar   PetscFunctionBegin;
4269566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
42733d57670SJed Brown   M          = A->rmap->N / bs;
428a6053eceSJunchao Zhang   mumps->val = aa->a;
429bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
4309371c9d4SSatish Balay     ai = aa->i;
4319371c9d4SSatish Balay     aj = aa->j;
43267877ebaSShri Abhyankar     nz = bs2 * aa->nz;
4339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
43467877ebaSShri Abhyankar     for (i = 0; i < M; i++) {
43567877ebaSShri Abhyankar       ajj = aj + ai[i];
43667877ebaSShri Abhyankar       rnz = ai[i + 1] - ai[i];
43767877ebaSShri Abhyankar       for (k = 0; k < rnz; k++) {
43867877ebaSShri Abhyankar         for (j = 0; j < bs; j++) {
43967877ebaSShri Abhyankar           for (m = 0; m < bs; m++) {
4409566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i * bs + m + shift, &row[idx]));
4419566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(bs * ajj[k] + j + shift, &col[idx]));
442a6053eceSJunchao Zhang             idx++;
44367877ebaSShri Abhyankar           }
44467877ebaSShri Abhyankar         }
44567877ebaSShri Abhyankar       }
44667877ebaSShri Abhyankar     }
447a6053eceSJunchao Zhang     mumps->irn = row;
448a6053eceSJunchao Zhang     mumps->jcn = col;
449a6053eceSJunchao Zhang     mumps->nnz = nz;
45067877ebaSShri Abhyankar   }
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45267877ebaSShri Abhyankar }
45367877ebaSShri Abhyankar 
454d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
455d71ae5a4SJacob Faibussowitsch {
45675480915SPierre Jolivet   const PetscInt *ai, *aj, *ajj;
457a6053eceSJunchao Zhang   PetscInt        bs;
458a6053eceSJunchao Zhang   PetscInt64      nz, rnz, i, j, k, m;
459a6053eceSJunchao Zhang   PetscMUMPSInt  *row, *col;
46075480915SPierre Jolivet   PetscScalar    *val;
46116ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa  = (Mat_SeqSBAIJ *)A->data;
46275480915SPierre Jolivet   const PetscInt  bs2 = aa->bs2, mbs = aa->mbs;
46338548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
464b94d7dedSBarry Smith   PetscBool isset, hermitian;
46538548759SBarry Smith #endif
46616ebf90aSShri Abhyankar 
46716ebf90aSShri Abhyankar   PetscFunctionBegin;
46838548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
469b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
470b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
47138548759SBarry Smith #endif
4722205254eSKarl Rupp   ai = aa->i;
4732205254eSKarl Rupp   aj = aa->j;
4749566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
47575480915SPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
476f3fa974cSJacob Faibussowitsch     const PetscInt64 alloc_size = aa->nz * bs2;
477f3fa974cSJacob Faibussowitsch 
478f3fa974cSJacob Faibussowitsch     PetscCall(PetscMalloc2(alloc_size, &row, alloc_size, &col));
479a6053eceSJunchao Zhang     if (bs > 1) {
480f3fa974cSJacob Faibussowitsch       PetscCall(PetscMalloc1(alloc_size, &mumps->val_alloc));
481a6053eceSJunchao Zhang       mumps->val = mumps->val_alloc;
48275480915SPierre Jolivet     } else {
483a6053eceSJunchao Zhang       mumps->val = aa->a;
48475480915SPierre Jolivet     }
485a6053eceSJunchao Zhang     mumps->irn = row;
486a6053eceSJunchao Zhang     mumps->jcn = col;
487a6053eceSJunchao Zhang   } else {
488a6053eceSJunchao Zhang     if (bs == 1) mumps->val = aa->a;
489a6053eceSJunchao Zhang     row = mumps->irn;
490a6053eceSJunchao Zhang     col = mumps->jcn;
491a6053eceSJunchao Zhang   }
492a6053eceSJunchao Zhang   val = mumps->val;
493185f6596SHong Zhang 
49416ebf90aSShri Abhyankar   nz = 0;
495a81fe166SPierre Jolivet   if (bs > 1) {
49675480915SPierre Jolivet     for (i = 0; i < mbs; i++) {
49716ebf90aSShri Abhyankar       rnz = ai[i + 1] - ai[i];
49867877ebaSShri Abhyankar       ajj = aj + ai[i];
49975480915SPierre Jolivet       for (j = 0; j < rnz; j++) {
50075480915SPierre Jolivet         for (k = 0; k < bs; k++) {
50175480915SPierre Jolivet           for (m = 0; m < bs; m++) {
502ec4f40fdSPierre Jolivet             if (ajj[j] > i || k >= m) {
50375480915SPierre Jolivet               if (reuse == MAT_INITIAL_MATRIX) {
5049566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(i * bs + m + shift, &row[nz]));
5059566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(ajj[j] * bs + k + shift, &col[nz]));
50675480915SPierre Jolivet               }
50775480915SPierre Jolivet               val[nz++] = aa->a[(ai[i] + j) * bs2 + m + k * bs];
50875480915SPierre Jolivet             }
50975480915SPierre Jolivet           }
51075480915SPierre Jolivet         }
51175480915SPierre Jolivet       }
51275480915SPierre Jolivet     }
513a81fe166SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) {
514a81fe166SPierre Jolivet     for (i = 0; i < mbs; i++) {
515a81fe166SPierre Jolivet       rnz = ai[i + 1] - ai[i];
516a81fe166SPierre Jolivet       ajj = aj + ai[i];
517a81fe166SPierre Jolivet       for (j = 0; j < rnz; j++) {
5189566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
5199566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[nz]));
520a6053eceSJunchao Zhang         nz++;
521a81fe166SPierre Jolivet       }
522a81fe166SPierre Jolivet     }
52308401ef6SPierre Jolivet     PetscCheck(nz == aa->nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Different numbers of nonzeros %" PetscInt64_FMT " != %" PetscInt_FMT, nz, aa->nz);
52475480915SPierre Jolivet   }
525a6053eceSJunchao Zhang   if (reuse == MAT_INITIAL_MATRIX) mumps->nnz = nz;
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52716ebf90aSShri Abhyankar }
52816ebf90aSShri Abhyankar 
529d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
530d71ae5a4SJacob Faibussowitsch {
53167877ebaSShri Abhyankar   const PetscInt    *ai, *aj, *ajj, *adiag, M = A->rmap->n;
532a6053eceSJunchao Zhang   PetscInt64         nz, rnz, i, j;
53367877ebaSShri Abhyankar   const PetscScalar *av, *v1;
53416ebf90aSShri Abhyankar   PetscScalar       *val;
535a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
536829b1710SHong Zhang   Mat_SeqAIJ        *aa = (Mat_SeqAIJ *)A->data;
53729b521d4Sstefano_zampini   PetscBool          missing;
53838548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
539b94d7dedSBarry Smith   PetscBool hermitian, isset;
54038548759SBarry Smith #endif
54116ebf90aSShri Abhyankar 
54216ebf90aSShri Abhyankar   PetscFunctionBegin;
54338548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
544b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
545b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
54638548759SBarry Smith #endif
5479566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(A, &av));
5489371c9d4SSatish Balay   ai    = aa->i;
5499371c9d4SSatish Balay   aj    = aa->j;
55016ebf90aSShri Abhyankar   adiag = aa->diag;
5519566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal_SeqAIJ(A, &missing, NULL));
552bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
5537ee00b23SStefano Zampini     /* count nz in the upper triangular part of A */
554829b1710SHong Zhang     nz = 0;
55529b521d4Sstefano_zampini     if (missing) {
55629b521d4Sstefano_zampini       for (i = 0; i < M; i++) {
55729b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i + 1])) {
55829b521d4Sstefano_zampini           for (j = ai[i]; j < ai[i + 1]; j++) {
55929b521d4Sstefano_zampini             if (aj[j] < i) continue;
56029b521d4Sstefano_zampini             nz++;
56129b521d4Sstefano_zampini           }
56229b521d4Sstefano_zampini         } else {
56329b521d4Sstefano_zampini           nz += ai[i + 1] - adiag[i];
56429b521d4Sstefano_zampini         }
56529b521d4Sstefano_zampini       }
56629b521d4Sstefano_zampini     } else {
567829b1710SHong Zhang       for (i = 0; i < M; i++) nz += ai[i + 1] - adiag[i];
56829b521d4Sstefano_zampini     }
5699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
5709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
571a6053eceSJunchao Zhang     mumps->nnz = nz;
572a6053eceSJunchao Zhang     mumps->irn = row;
573a6053eceSJunchao Zhang     mumps->jcn = col;
574a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
575185f6596SHong Zhang 
57616ebf90aSShri Abhyankar     nz = 0;
57729b521d4Sstefano_zampini     if (missing) {
57829b521d4Sstefano_zampini       for (i = 0; i < M; i++) {
57929b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i + 1])) {
58029b521d4Sstefano_zampini           for (j = ai[i]; j < ai[i + 1]; j++) {
58129b521d4Sstefano_zampini             if (aj[j] < i) continue;
5829566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
5839566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(aj[j] + shift, &col[nz]));
58429b521d4Sstefano_zampini             val[nz] = av[j];
58529b521d4Sstefano_zampini             nz++;
58629b521d4Sstefano_zampini           }
58729b521d4Sstefano_zampini         } else {
58829b521d4Sstefano_zampini           rnz = ai[i + 1] - adiag[i];
58929b521d4Sstefano_zampini           ajj = aj + adiag[i];
59029b521d4Sstefano_zampini           v1  = av + adiag[i];
59129b521d4Sstefano_zampini           for (j = 0; j < rnz; j++) {
5929566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
5939566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[nz]));
594a6053eceSJunchao Zhang             val[nz++] = v1[j];
59529b521d4Sstefano_zampini           }
59629b521d4Sstefano_zampini         }
59729b521d4Sstefano_zampini       }
59829b521d4Sstefano_zampini     } else {
59916ebf90aSShri Abhyankar       for (i = 0; i < M; i++) {
60016ebf90aSShri Abhyankar         rnz = ai[i + 1] - adiag[i];
60167877ebaSShri Abhyankar         ajj = aj + adiag[i];
602cf3759fdSShri Abhyankar         v1  = av + adiag[i];
60367877ebaSShri Abhyankar         for (j = 0; j < rnz; j++) {
6049566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
6059566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[nz]));
606a6053eceSJunchao Zhang           val[nz++] = v1[j];
60716ebf90aSShri Abhyankar         }
60816ebf90aSShri Abhyankar       }
60929b521d4Sstefano_zampini     }
610397b6df1SKris Buschelman   } else {
611a6053eceSJunchao Zhang     nz  = 0;
612a6053eceSJunchao Zhang     val = mumps->val;
61329b521d4Sstefano_zampini     if (missing) {
61416ebf90aSShri Abhyankar       for (i = 0; i < M; i++) {
61529b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i + 1])) {
61629b521d4Sstefano_zampini           for (j = ai[i]; j < ai[i + 1]; j++) {
61729b521d4Sstefano_zampini             if (aj[j] < i) continue;
61829b521d4Sstefano_zampini             val[nz++] = av[j];
61929b521d4Sstefano_zampini           }
62029b521d4Sstefano_zampini         } else {
62116ebf90aSShri Abhyankar           rnz = ai[i + 1] - adiag[i];
62267877ebaSShri Abhyankar           v1  = av + adiag[i];
623ad540459SPierre Jolivet           for (j = 0; j < rnz; j++) val[nz++] = v1[j];
62416ebf90aSShri Abhyankar         }
62516ebf90aSShri Abhyankar       }
62629b521d4Sstefano_zampini     } else {
62716ebf90aSShri Abhyankar       for (i = 0; i < M; i++) {
62816ebf90aSShri Abhyankar         rnz = ai[i + 1] - adiag[i];
62916ebf90aSShri Abhyankar         v1  = av + adiag[i];
630ad540459SPierre Jolivet         for (j = 0; j < rnz; j++) val[nz++] = v1[j];
63116ebf90aSShri Abhyankar       }
63216ebf90aSShri Abhyankar     }
63329b521d4Sstefano_zampini   }
6349566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63616ebf90aSShri Abhyankar }
63716ebf90aSShri Abhyankar 
638d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
639d71ae5a4SJacob Faibussowitsch {
640a6053eceSJunchao Zhang   const PetscInt    *ai, *aj, *bi, *bj, *garray, *ajj, *bjj;
641a6053eceSJunchao Zhang   PetscInt           bs;
642a6053eceSJunchao Zhang   PetscInt64         rstart, nz, i, j, k, m, jj, irow, countA, countB;
643a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
64416ebf90aSShri Abhyankar   const PetscScalar *av, *bv, *v1, *v2;
64516ebf90aSShri Abhyankar   PetscScalar       *val;
646397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ *)A->data;
647397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ *)(mat->A)->data;
648397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ *)(mat->B)->data;
649ec4f40fdSPierre Jolivet   const PetscInt     bs2 = aa->bs2, mbs = aa->mbs;
65038548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
651b94d7dedSBarry Smith   PetscBool hermitian, isset;
65238548759SBarry Smith #endif
65316ebf90aSShri Abhyankar 
65416ebf90aSShri Abhyankar   PetscFunctionBegin;
65538548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
656b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
657b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
65838548759SBarry Smith #endif
6599566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
66038548759SBarry Smith   rstart = A->rmap->rstart;
66138548759SBarry Smith   ai     = aa->i;
66238548759SBarry Smith   aj     = aa->j;
66338548759SBarry Smith   bi     = bb->i;
66438548759SBarry Smith   bj     = bb->j;
66538548759SBarry Smith   av     = aa->a;
66638548759SBarry Smith   bv     = bb->a;
667397b6df1SKris Buschelman 
6682205254eSKarl Rupp   garray = mat->garray;
6692205254eSKarl Rupp 
670bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
671a6053eceSJunchao Zhang     nz = (aa->nz + bb->nz) * bs2; /* just a conservative estimate */
6729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
6739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
674a6053eceSJunchao Zhang     /* can not decide the exact mumps->nnz now because of the SBAIJ */
675a6053eceSJunchao Zhang     mumps->irn = row;
676a6053eceSJunchao Zhang     mumps->jcn = col;
677a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
678397b6df1SKris Buschelman   } else {
679a6053eceSJunchao Zhang     val = mumps->val;
680397b6df1SKris Buschelman   }
681397b6df1SKris Buschelman 
6829371c9d4SSatish Balay   jj   = 0;
6839371c9d4SSatish Balay   irow = rstart;
684ec4f40fdSPierre Jolivet   for (i = 0; i < mbs; i++) {
685397b6df1SKris Buschelman     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
686397b6df1SKris Buschelman     countA = ai[i + 1] - ai[i];
687397b6df1SKris Buschelman     countB = bi[i + 1] - bi[i];
688397b6df1SKris Buschelman     bjj    = bj + bi[i];
689ec4f40fdSPierre Jolivet     v1     = av + ai[i] * bs2;
690ec4f40fdSPierre Jolivet     v2     = bv + bi[i] * bs2;
691397b6df1SKris Buschelman 
692ec4f40fdSPierre Jolivet     if (bs > 1) {
693ec4f40fdSPierre Jolivet       /* A-part */
694ec4f40fdSPierre Jolivet       for (j = 0; j < countA; j++) {
695ec4f40fdSPierre Jolivet         for (k = 0; k < bs; k++) {
696ec4f40fdSPierre Jolivet           for (m = 0; m < bs; m++) {
697ec4f40fdSPierre Jolivet             if (rstart + ajj[j] * bs > irow || k >= m) {
698ec4f40fdSPierre Jolivet               if (reuse == MAT_INITIAL_MATRIX) {
6999566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(irow + m + shift, &row[jj]));
7009566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(rstart + ajj[j] * bs + k + shift, &col[jj]));
701ec4f40fdSPierre Jolivet               }
702ec4f40fdSPierre Jolivet               val[jj++] = v1[j * bs2 + m + k * bs];
703ec4f40fdSPierre Jolivet             }
704ec4f40fdSPierre Jolivet           }
705ec4f40fdSPierre Jolivet         }
706ec4f40fdSPierre Jolivet       }
707ec4f40fdSPierre Jolivet 
708ec4f40fdSPierre Jolivet       /* B-part */
709ec4f40fdSPierre Jolivet       for (j = 0; j < countB; j++) {
710ec4f40fdSPierre Jolivet         for (k = 0; k < bs; k++) {
711ec4f40fdSPierre Jolivet           for (m = 0; m < bs; m++) {
712ec4f40fdSPierre Jolivet             if (reuse == MAT_INITIAL_MATRIX) {
7139566063dSJacob Faibussowitsch               PetscCall(PetscMUMPSIntCast(irow + m + shift, &row[jj]));
7149566063dSJacob Faibussowitsch               PetscCall(PetscMUMPSIntCast(garray[bjj[j]] * bs + k + shift, &col[jj]));
715ec4f40fdSPierre Jolivet             }
716ec4f40fdSPierre Jolivet             val[jj++] = v2[j * bs2 + m + k * bs];
717ec4f40fdSPierre Jolivet           }
718ec4f40fdSPierre Jolivet         }
719ec4f40fdSPierre Jolivet       }
720ec4f40fdSPierre Jolivet     } else {
721397b6df1SKris Buschelman       /* A-part */
722397b6df1SKris Buschelman       for (j = 0; j < countA; j++) {
723bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
7249566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
7259566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift, &col[jj]));
726397b6df1SKris Buschelman         }
72716ebf90aSShri Abhyankar         val[jj++] = v1[j];
728397b6df1SKris Buschelman       }
72916ebf90aSShri Abhyankar 
73016ebf90aSShri Abhyankar       /* B-part */
73116ebf90aSShri Abhyankar       for (j = 0; j < countB; j++) {
732bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
7339566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
7349566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift, &col[jj]));
735397b6df1SKris Buschelman         }
73616ebf90aSShri Abhyankar         val[jj++] = v2[j];
73716ebf90aSShri Abhyankar       }
73816ebf90aSShri Abhyankar     }
739ec4f40fdSPierre Jolivet     irow += bs;
740ec4f40fdSPierre Jolivet   }
741a6053eceSJunchao Zhang   mumps->nnz = jj;
7423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74316ebf90aSShri Abhyankar }
74416ebf90aSShri Abhyankar 
745d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
746d71ae5a4SJacob Faibussowitsch {
74716ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj, *garray, m = A->rmap->n, *ajj, *bjj;
748a6053eceSJunchao Zhang   PetscInt64         rstart, nz, i, j, jj, irow, countA, countB;
749a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
75016ebf90aSShri Abhyankar   const PetscScalar *av, *bv, *v1, *v2;
75116ebf90aSShri Abhyankar   PetscScalar       *val;
752a3d589ffSStefano Zampini   Mat                Ad, Ao;
753a3d589ffSStefano Zampini   Mat_SeqAIJ        *aa;
754a3d589ffSStefano Zampini   Mat_SeqAIJ        *bb;
75516ebf90aSShri Abhyankar 
75616ebf90aSShri Abhyankar   PetscFunctionBegin;
7579566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A, &Ad, &Ao, &garray));
7589566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ad, &av));
7599566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ao, &bv));
760a3d589ffSStefano Zampini 
761a3d589ffSStefano Zampini   aa = (Mat_SeqAIJ *)(Ad)->data;
762a3d589ffSStefano Zampini   bb = (Mat_SeqAIJ *)(Ao)->data;
76338548759SBarry Smith   ai = aa->i;
76438548759SBarry Smith   aj = aa->j;
76538548759SBarry Smith   bi = bb->i;
76638548759SBarry Smith   bj = bb->j;
76716ebf90aSShri Abhyankar 
768a3d589ffSStefano Zampini   rstart = A->rmap->rstart;
7692205254eSKarl Rupp 
770bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
771a6053eceSJunchao Zhang     nz = (PetscInt64)aa->nz + bb->nz; /* make sure the sum won't overflow PetscInt */
7729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
7739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
774a6053eceSJunchao Zhang     mumps->nnz = nz;
775a6053eceSJunchao Zhang     mumps->irn = row;
776a6053eceSJunchao Zhang     mumps->jcn = col;
777a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
77816ebf90aSShri Abhyankar   } else {
779a6053eceSJunchao Zhang     val = mumps->val;
78016ebf90aSShri Abhyankar   }
78116ebf90aSShri Abhyankar 
7829371c9d4SSatish Balay   jj   = 0;
7839371c9d4SSatish Balay   irow = rstart;
78416ebf90aSShri Abhyankar   for (i = 0; i < m; i++) {
78516ebf90aSShri Abhyankar     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
78616ebf90aSShri Abhyankar     countA = ai[i + 1] - ai[i];
78716ebf90aSShri Abhyankar     countB = bi[i + 1] - bi[i];
78816ebf90aSShri Abhyankar     bjj    = bj + bi[i];
78916ebf90aSShri Abhyankar     v1     = av + ai[i];
79016ebf90aSShri Abhyankar     v2     = bv + bi[i];
79116ebf90aSShri Abhyankar 
79216ebf90aSShri Abhyankar     /* A-part */
79316ebf90aSShri Abhyankar     for (j = 0; j < countA; j++) {
794bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
7959566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
7969566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift, &col[jj]));
79716ebf90aSShri Abhyankar       }
79816ebf90aSShri Abhyankar       val[jj++] = v1[j];
79916ebf90aSShri Abhyankar     }
80016ebf90aSShri Abhyankar 
80116ebf90aSShri Abhyankar     /* B-part */
80216ebf90aSShri Abhyankar     for (j = 0; j < countB; j++) {
803bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
8049566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
8059566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift, &col[jj]));
80616ebf90aSShri Abhyankar       }
80716ebf90aSShri Abhyankar       val[jj++] = v2[j];
80816ebf90aSShri Abhyankar     }
80916ebf90aSShri Abhyankar     irow++;
81016ebf90aSShri Abhyankar   }
8119566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ad, &av));
8129566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ao, &bv));
8133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81416ebf90aSShri Abhyankar }
81516ebf90aSShri Abhyankar 
816d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
817d71ae5a4SJacob Faibussowitsch {
81867877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat = (Mat_MPIBAIJ *)A->data;
81967877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa  = (Mat_SeqBAIJ *)(mat->A)->data;
82067877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ *)(mat->B)->data;
82167877ebaSShri Abhyankar   const PetscInt    *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j, *ajj, *bjj;
822d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray, mbs = mat->mbs, rstart = A->rmap->rstart;
82333d57670SJed Brown   const PetscInt     bs2 = mat->bs2;
824a6053eceSJunchao Zhang   PetscInt           bs;
825a6053eceSJunchao Zhang   PetscInt64         nz, i, j, k, n, jj, irow, countA, countB, idx;
826a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
82767877ebaSShri Abhyankar   const PetscScalar *av = aa->a, *bv = bb->a, *v1, *v2;
82867877ebaSShri Abhyankar   PetscScalar       *val;
82967877ebaSShri Abhyankar 
83067877ebaSShri Abhyankar   PetscFunctionBegin;
8319566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
832bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
83367877ebaSShri Abhyankar     nz = bs2 * (aa->nz + bb->nz);
8349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
8359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
836a6053eceSJunchao Zhang     mumps->nnz = nz;
837a6053eceSJunchao Zhang     mumps->irn = row;
838a6053eceSJunchao Zhang     mumps->jcn = col;
839a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
84067877ebaSShri Abhyankar   } else {
841a6053eceSJunchao Zhang     val = mumps->val;
84267877ebaSShri Abhyankar   }
84367877ebaSShri Abhyankar 
8449371c9d4SSatish Balay   jj   = 0;
8459371c9d4SSatish Balay   irow = rstart;
84667877ebaSShri Abhyankar   for (i = 0; i < mbs; i++) {
84767877ebaSShri Abhyankar     countA = ai[i + 1] - ai[i];
84867877ebaSShri Abhyankar     countB = bi[i + 1] - bi[i];
84967877ebaSShri Abhyankar     ajj    = aj + ai[i];
85067877ebaSShri Abhyankar     bjj    = bj + bi[i];
85167877ebaSShri Abhyankar     v1     = av + bs2 * ai[i];
85267877ebaSShri Abhyankar     v2     = bv + bs2 * bi[i];
85367877ebaSShri Abhyankar 
85467877ebaSShri Abhyankar     idx = 0;
85567877ebaSShri Abhyankar     /* A-part */
85667877ebaSShri Abhyankar     for (k = 0; k < countA; k++) {
85767877ebaSShri Abhyankar       for (j = 0; j < bs; j++) {
85867877ebaSShri Abhyankar         for (n = 0; n < bs; n++) {
859bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
8609566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(irow + n + shift, &row[jj]));
8619566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(rstart + bs * ajj[k] + j + shift, &col[jj]));
86267877ebaSShri Abhyankar           }
86367877ebaSShri Abhyankar           val[jj++] = v1[idx++];
86467877ebaSShri Abhyankar         }
86567877ebaSShri Abhyankar       }
86667877ebaSShri Abhyankar     }
86767877ebaSShri Abhyankar 
86867877ebaSShri Abhyankar     idx = 0;
86967877ebaSShri Abhyankar     /* B-part */
87067877ebaSShri Abhyankar     for (k = 0; k < countB; k++) {
87167877ebaSShri Abhyankar       for (j = 0; j < bs; j++) {
87267877ebaSShri Abhyankar         for (n = 0; n < bs; n++) {
873bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
8749566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(irow + n + shift, &row[jj]));
8759566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(bs * garray[bjj[k]] + j + shift, &col[jj]));
87667877ebaSShri Abhyankar           }
877d985c460SShri Abhyankar           val[jj++] = v2[idx++];
87867877ebaSShri Abhyankar         }
87967877ebaSShri Abhyankar       }
88067877ebaSShri Abhyankar     }
881d985c460SShri Abhyankar     irow += bs;
88267877ebaSShri Abhyankar   }
8833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88467877ebaSShri Abhyankar }
88567877ebaSShri Abhyankar 
886d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
887d71ae5a4SJacob Faibussowitsch {
88816ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *adiag, *bi, *bj, *garray, m = A->rmap->n, *ajj, *bjj;
889a6053eceSJunchao Zhang   PetscInt64         rstart, nz, nza, nzb, i, j, jj, irow, countA, countB;
890a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
89116ebf90aSShri Abhyankar   const PetscScalar *av, *bv, *v1, *v2;
89216ebf90aSShri Abhyankar   PetscScalar       *val;
893a3d589ffSStefano Zampini   Mat                Ad, Ao;
894a3d589ffSStefano Zampini   Mat_SeqAIJ        *aa;
895a3d589ffSStefano Zampini   Mat_SeqAIJ        *bb;
89638548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
897b94d7dedSBarry Smith   PetscBool hermitian, isset;
89838548759SBarry Smith #endif
89916ebf90aSShri Abhyankar 
90016ebf90aSShri Abhyankar   PetscFunctionBegin;
90138548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
902b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
903b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
90438548759SBarry Smith #endif
9059566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A, &Ad, &Ao, &garray));
9069566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ad, &av));
9079566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ao, &bv));
908a3d589ffSStefano Zampini 
909a3d589ffSStefano Zampini   aa    = (Mat_SeqAIJ *)(Ad)->data;
910a3d589ffSStefano Zampini   bb    = (Mat_SeqAIJ *)(Ao)->data;
91138548759SBarry Smith   ai    = aa->i;
91238548759SBarry Smith   aj    = aa->j;
91338548759SBarry Smith   adiag = aa->diag;
91438548759SBarry Smith   bi    = bb->i;
91538548759SBarry Smith   bj    = bb->j;
9162205254eSKarl Rupp 
91716ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
91816ebf90aSShri Abhyankar 
919bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
920e0bace9bSHong Zhang     nza = 0; /* num of upper triangular entries in mat->A, including diagonals */
921e0bace9bSHong Zhang     nzb = 0; /* num of upper triangular entries in mat->B */
92216ebf90aSShri Abhyankar     for (i = 0; i < m; i++) {
923e0bace9bSHong Zhang       nza += (ai[i + 1] - adiag[i]);
92416ebf90aSShri Abhyankar       countB = bi[i + 1] - bi[i];
92516ebf90aSShri Abhyankar       bjj    = bj + bi[i];
926e0bace9bSHong Zhang       for (j = 0; j < countB; j++) {
927e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
928e0bace9bSHong Zhang       }
929e0bace9bSHong Zhang     }
93016ebf90aSShri Abhyankar 
931e0bace9bSHong Zhang     nz = nza + nzb; /* total nz of upper triangular part of mat */
9329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
9339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
934a6053eceSJunchao Zhang     mumps->nnz = nz;
935a6053eceSJunchao Zhang     mumps->irn = row;
936a6053eceSJunchao Zhang     mumps->jcn = col;
937a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
93816ebf90aSShri Abhyankar   } else {
939a6053eceSJunchao Zhang     val = mumps->val;
94016ebf90aSShri Abhyankar   }
94116ebf90aSShri Abhyankar 
9429371c9d4SSatish Balay   jj   = 0;
9439371c9d4SSatish Balay   irow = rstart;
94416ebf90aSShri Abhyankar   for (i = 0; i < m; i++) {
94516ebf90aSShri Abhyankar     ajj    = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */
94616ebf90aSShri Abhyankar     v1     = av + adiag[i];
94716ebf90aSShri Abhyankar     countA = ai[i + 1] - adiag[i];
94816ebf90aSShri Abhyankar     countB = bi[i + 1] - bi[i];
94916ebf90aSShri Abhyankar     bjj    = bj + bi[i];
95016ebf90aSShri Abhyankar     v2     = bv + bi[i];
95116ebf90aSShri Abhyankar 
95216ebf90aSShri Abhyankar     /* A-part */
95316ebf90aSShri Abhyankar     for (j = 0; j < countA; j++) {
954bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
9559566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
9569566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift, &col[jj]));
95716ebf90aSShri Abhyankar       }
95816ebf90aSShri Abhyankar       val[jj++] = v1[j];
95916ebf90aSShri Abhyankar     }
96016ebf90aSShri Abhyankar 
96116ebf90aSShri Abhyankar     /* B-part */
96216ebf90aSShri Abhyankar     for (j = 0; j < countB; j++) {
96316ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
964bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
9659566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
9669566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift, &col[jj]));
96716ebf90aSShri Abhyankar         }
96816ebf90aSShri Abhyankar         val[jj++] = v2[j];
96916ebf90aSShri Abhyankar       }
970397b6df1SKris Buschelman     }
971397b6df1SKris Buschelman     irow++;
972397b6df1SKris Buschelman   }
9739566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ad, &av));
9749566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ao, &bv));
9753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
976397b6df1SKris Buschelman }
977397b6df1SKris Buschelman 
978d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_MUMPS(Mat A)
979d71ae5a4SJacob Faibussowitsch {
980a6053eceSJunchao Zhang   Mat_MUMPS *mumps = (Mat_MUMPS *)A->data;
981b24902e0SBarry Smith 
982397b6df1SKris Buschelman   PetscFunctionBegin;
9839566063dSJacob Faibussowitsch   PetscCall(PetscFree2(mumps->id.sol_loc, mumps->id.isol_loc));
9849566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mumps->scat_rhs));
9859566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mumps->scat_sol));
9869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mumps->b_seq));
9879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mumps->x_seq));
9889566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.perm_in));
9899566063dSJacob Faibussowitsch   PetscCall(PetscFree2(mumps->irn, mumps->jcn));
9909566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->val_alloc));
9919566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->info));
992413bcc21SPierre Jolivet   PetscCall(PetscFree(mumps->ICNTL_pre));
993413bcc21SPierre Jolivet   PetscCall(PetscFree(mumps->CNTL_pre));
9949566063dSJacob Faibussowitsch   PetscCall(MatMumpsResetSchur_Private(mumps));
995413bcc21SPierre Jolivet   if (mumps->id.job != JOB_NULL) { /* cannot call PetscMUMPS_c() if JOB_INIT has never been called for this instance */
996a5e57a09SHong Zhang     mumps->id.job = JOB_END;
9973ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
99808401ef6SPierre Jolivet     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in MatDestroy_MUMPS: INFOG(1)=%d", mumps->id.INFOG(1));
999413bcc21SPierre Jolivet     if (mumps->mumps_comm != MPI_COMM_NULL) {
1000413bcc21SPierre Jolivet       if (PetscDefined(HAVE_OPENMP_SUPPORT) && mumps->use_petsc_omp_support) PetscCallMPI(MPI_Comm_free(&mumps->mumps_comm));
1001413bcc21SPierre Jolivet       else PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)A), &mumps->mumps_comm));
1002413bcc21SPierre Jolivet     }
1003413bcc21SPierre Jolivet   }
10043ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
100567602552SJunchao Zhang   if (mumps->use_petsc_omp_support) {
10069566063dSJacob Faibussowitsch     PetscCall(PetscOmpCtrlDestroy(&mumps->omp_ctrl));
10079566063dSJacob Faibussowitsch     PetscCall(PetscFree2(mumps->rhs_loc, mumps->rhs_recvbuf));
10089566063dSJacob Faibussowitsch     PetscCall(PetscFree3(mumps->rhs_nrow, mumps->rhs_recvcounts, mumps->rhs_disps));
100967602552SJunchao Zhang   }
10103ab56b82SJunchao Zhang #endif
10119566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->ia_alloc));
10129566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->ja_alloc));
10139566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->recvcount));
10149566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->reqs));
10159566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->irhs_loc));
10169566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1017bf0cc555SLisandro Dalcin 
101897969023SHong Zhang   /* clear composed functions */
10199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
10209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorSetSchurIS_C", NULL));
10219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorCreateSchurComplement_C", NULL));
10229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsSetIcntl_C", NULL));
10239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetIcntl_C", NULL));
10249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsSetCntl_C", NULL));
10259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetCntl_C", NULL));
10269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInfo_C", NULL));
10279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInfog_C", NULL));
10289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetRinfo_C", NULL));
10299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetRinfog_C", NULL));
10305c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetNullPivots_C", NULL));
10319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInverse_C", NULL));
10329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInverseTranspose_C", NULL));
10333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1034397b6df1SKris Buschelman }
1035397b6df1SKris Buschelman 
103667602552SJunchao Zhang /* Set up the distributed RHS info for MUMPS. <nrhs> is the number of RHS. <array> points to start of RHS on the local processor. */
1037d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsSetUpDistRHSInfo(Mat A, PetscInt nrhs, const PetscScalar *array)
1038d71ae5a4SJacob Faibussowitsch {
103967602552SJunchao Zhang   Mat_MUMPS        *mumps   = (Mat_MUMPS *)A->data;
104067602552SJunchao Zhang   const PetscMPIInt ompsize = mumps->omp_comm_size;
104167602552SJunchao Zhang   PetscInt          i, m, M, rstart;
104267602552SJunchao Zhang 
104367602552SJunchao Zhang   PetscFunctionBegin;
10449566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, NULL));
10459566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
104608401ef6SPierre Jolivet   PetscCheck(M <= PETSC_MUMPS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "PetscInt too long for PetscMUMPSInt");
104767602552SJunchao Zhang   if (ompsize == 1) {
104867602552SJunchao Zhang     if (!mumps->irhs_loc) {
104967602552SJunchao Zhang       mumps->nloc_rhs = m;
10509566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &mumps->irhs_loc));
10519566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
105267602552SJunchao Zhang       for (i = 0; i < m; i++) mumps->irhs_loc[i] = rstart + i + 1; /* use 1-based indices */
105367602552SJunchao Zhang     }
105467602552SJunchao Zhang     mumps->id.rhs_loc = (MumpsScalar *)array;
105567602552SJunchao Zhang   } else {
105667602552SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
105767602552SJunchao Zhang     const PetscInt *ranges;
105867602552SJunchao Zhang     PetscMPIInt     j, k, sendcount, *petsc_ranks, *omp_ranks;
105967602552SJunchao Zhang     MPI_Group       petsc_group, omp_group;
106067602552SJunchao Zhang     PetscScalar    *recvbuf = NULL;
106167602552SJunchao Zhang 
106267602552SJunchao Zhang     if (mumps->is_omp_master) {
106367602552SJunchao Zhang       /* Lazily initialize the omp stuff for distributed rhs */
106467602552SJunchao Zhang       if (!mumps->irhs_loc) {
10659566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(ompsize, &omp_ranks, ompsize, &petsc_ranks));
10669566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(ompsize, &mumps->rhs_nrow, ompsize, &mumps->rhs_recvcounts, ompsize, &mumps->rhs_disps));
10679566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_group(mumps->petsc_comm, &petsc_group));
10689566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_group(mumps->omp_comm, &omp_group));
106967602552SJunchao Zhang         for (j = 0; j < ompsize; j++) omp_ranks[j] = j;
10709566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_translate_ranks(omp_group, ompsize, omp_ranks, petsc_group, petsc_ranks));
107167602552SJunchao Zhang 
107267602552SJunchao Zhang         /* Populate mumps->irhs_loc[], rhs_nrow[] */
107367602552SJunchao Zhang         mumps->nloc_rhs = 0;
10749566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRanges(A, &ranges));
107567602552SJunchao Zhang         for (j = 0; j < ompsize; j++) {
107667602552SJunchao Zhang           mumps->rhs_nrow[j] = ranges[petsc_ranks[j] + 1] - ranges[petsc_ranks[j]];
107767602552SJunchao Zhang           mumps->nloc_rhs += mumps->rhs_nrow[j];
107867602552SJunchao Zhang         }
10799566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(mumps->nloc_rhs, &mumps->irhs_loc));
108067602552SJunchao Zhang         for (j = k = 0; j < ompsize; j++) {
108167602552SJunchao Zhang           for (i = ranges[petsc_ranks[j]]; i < ranges[petsc_ranks[j] + 1]; i++, k++) mumps->irhs_loc[k] = i + 1; /* uses 1-based indices */
108267602552SJunchao Zhang         }
108367602552SJunchao Zhang 
10849566063dSJacob Faibussowitsch         PetscCall(PetscFree2(omp_ranks, petsc_ranks));
10859566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_free(&petsc_group));
10869566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_free(&omp_group));
108767602552SJunchao Zhang       }
108867602552SJunchao Zhang 
108967602552SJunchao Zhang       /* Realloc buffers when current nrhs is bigger than what we have met */
109067602552SJunchao Zhang       if (nrhs > mumps->max_nrhs) {
10919566063dSJacob Faibussowitsch         PetscCall(PetscFree2(mumps->rhs_loc, mumps->rhs_recvbuf));
10929566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(mumps->nloc_rhs * nrhs, &mumps->rhs_loc, mumps->nloc_rhs * nrhs, &mumps->rhs_recvbuf));
109367602552SJunchao Zhang         mumps->max_nrhs = nrhs;
109467602552SJunchao Zhang       }
109567602552SJunchao Zhang 
109667602552SJunchao Zhang       /* Setup recvcounts[], disps[], recvbuf on omp rank 0 for the upcoming MPI_Gatherv */
10979566063dSJacob Faibussowitsch       for (j = 0; j < ompsize; j++) PetscCall(PetscMPIIntCast(mumps->rhs_nrow[j] * nrhs, &mumps->rhs_recvcounts[j]));
109867602552SJunchao Zhang       mumps->rhs_disps[0] = 0;
109967602552SJunchao Zhang       for (j = 1; j < ompsize; j++) {
110067602552SJunchao Zhang         mumps->rhs_disps[j] = mumps->rhs_disps[j - 1] + mumps->rhs_recvcounts[j - 1];
110108401ef6SPierre Jolivet         PetscCheck(mumps->rhs_disps[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "PetscMPIInt overflow!");
110267602552SJunchao Zhang       }
110367602552SJunchao Zhang       recvbuf = (nrhs == 1) ? mumps->rhs_loc : mumps->rhs_recvbuf; /* Directly use rhs_loc[] as recvbuf. Single rhs is common in Ax=b */
110467602552SJunchao Zhang     }
110567602552SJunchao Zhang 
11069566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(m * nrhs, &sendcount));
11079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(array, sendcount, MPIU_SCALAR, recvbuf, mumps->rhs_recvcounts, mumps->rhs_disps, MPIU_SCALAR, 0, mumps->omp_comm));
110867602552SJunchao Zhang 
110967602552SJunchao Zhang     if (mumps->is_omp_master) {
111067602552SJunchao Zhang       if (nrhs > 1) { /* Copy & re-arrange data from rhs_recvbuf[] to mumps->rhs_loc[] only when there are multiple rhs */
111167602552SJunchao Zhang         PetscScalar *dst, *dstbase = mumps->rhs_loc;
111267602552SJunchao Zhang         for (j = 0; j < ompsize; j++) {
111367602552SJunchao Zhang           const PetscScalar *src = mumps->rhs_recvbuf + mumps->rhs_disps[j];
111467602552SJunchao Zhang           dst                    = dstbase;
111567602552SJunchao Zhang           for (i = 0; i < nrhs; i++) {
11169566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(dst, src, mumps->rhs_nrow[j]));
111767602552SJunchao Zhang             src += mumps->rhs_nrow[j];
111867602552SJunchao Zhang             dst += mumps->nloc_rhs;
111967602552SJunchao Zhang           }
112067602552SJunchao Zhang           dstbase += mumps->rhs_nrow[j];
112167602552SJunchao Zhang         }
112267602552SJunchao Zhang       }
112367602552SJunchao Zhang       mumps->id.rhs_loc = (MumpsScalar *)mumps->rhs_loc;
112467602552SJunchao Zhang     }
112567602552SJunchao Zhang #endif /* PETSC_HAVE_OPENMP_SUPPORT */
112667602552SJunchao Zhang   }
112767602552SJunchao Zhang   mumps->id.nrhs     = nrhs;
112867602552SJunchao Zhang   mumps->id.nloc_rhs = mumps->nloc_rhs;
112967602552SJunchao Zhang   mumps->id.lrhs_loc = mumps->nloc_rhs;
113067602552SJunchao Zhang   mumps->id.irhs_loc = mumps->irhs_loc;
11313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113267602552SJunchao Zhang }
113367602552SJunchao Zhang 
1134d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_MUMPS(Mat A, Vec b, Vec x)
1135d71ae5a4SJacob Faibussowitsch {
1136e69c285eSBarry Smith   Mat_MUMPS         *mumps  = (Mat_MUMPS *)A->data;
113725aac85cSJunchao Zhang   const PetscScalar *rarray = NULL;
1138d54de34fSKris Buschelman   PetscScalar       *array;
1139329ec9b3SHong Zhang   IS                 is_iden, is_petsc;
1140329ec9b3SHong Zhang   PetscInt           i;
1141cc86f929SStefano Zampini   PetscBool          second_solve = PETSC_FALSE;
1142883f2eb9SBarry Smith   static PetscBool   cite1 = PETSC_FALSE, cite2 = PETSC_FALSE;
1143397b6df1SKris Buschelman 
1144397b6df1SKris Buschelman   PetscFunctionBegin;
11459371c9d4SSatish Balay   PetscCall(PetscCitationsRegister("@article{MUMPS01,\n  author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n  title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n  journal = {SIAM "
11469371c9d4SSatish Balay                                    "Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",
11479371c9d4SSatish Balay                                    &cite1));
11489371c9d4SSatish Balay   PetscCall(PetscCitationsRegister("@article{MUMPS02,\n  author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n  title = {Hybrid scheduling for the parallel solution of linear systems},\n  journal = {Parallel "
11499371c9d4SSatish Balay                                    "Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",
11509371c9d4SSatish Balay                                    &cite2));
11512aca8efcSHong Zhang 
1152603e8f96SBarry Smith   if (A->factorerrortype) {
11539566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
11549566063dSJacob Faibussowitsch     PetscCall(VecSetInf(x));
11553ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
11562aca8efcSHong Zhang   }
11572aca8efcSHong Zhang 
1158a5e57a09SHong Zhang   mumps->id.nrhs = 1;
11592d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
116025aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
116167602552SJunchao Zhang       mumps->id.ICNTL(20) = 10; /* dense distributed RHS */
11629566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(b, &rarray));
11639566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetUpDistRHSInfo(A, 1, rarray));
116425aac85cSJunchao Zhang     } else {
116541ffd417SStefano Zampini       mumps->id.ICNTL(20) = 0; /* dense centralized RHS; Scatter b into a sequential rhs vector*/
11669566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mumps->scat_rhs, b, mumps->b_seq, INSERT_VALUES, SCATTER_FORWARD));
11679566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mumps->scat_rhs, b, mumps->b_seq, INSERT_VALUES, SCATTER_FORWARD));
116867602552SJunchao Zhang       if (!mumps->myid) {
11699566063dSJacob Faibussowitsch         PetscCall(VecGetArray(mumps->b_seq, &array));
117067602552SJunchao Zhang         mumps->id.rhs = (MumpsScalar *)array;
117167602552SJunchao Zhang       }
117225aac85cSJunchao Zhang     }
11733ab56b82SJunchao Zhang   } else {                   /* petsc_size == 1 */
117467602552SJunchao Zhang     mumps->id.ICNTL(20) = 0; /* dense centralized RHS */
11759566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, x));
11769566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &array));
1177940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar *)array;
1178397b6df1SKris Buschelman   }
1179397b6df1SKris Buschelman 
1180cc86f929SStefano Zampini   /*
1181cc86f929SStefano Zampini      handle condensation step of Schur complement (if any)
1182cc86f929SStefano Zampini      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
1183cc86f929SStefano Zampini      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
1184cc86f929SStefano Zampini      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
1185cc86f929SStefano Zampini      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
1186cc86f929SStefano Zampini   */
1187583f777eSStefano Zampini   if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
118808401ef6SPierre Jolivet     PetscCheck(mumps->petsc_size <= 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Parallel Schur complements not yet supported from PETSc");
1189cc86f929SStefano Zampini     second_solve = PETSC_TRUE;
11909566063dSJacob Faibussowitsch     PetscCall(MatMumpsHandleSchur_Private(A, PETSC_FALSE));
1191cc86f929SStefano Zampini   }
1192397b6df1SKris Buschelman   /* solve phase */
1193a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
11943ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
119508401ef6SPierre Jolivet   PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in solve phase: INFOG(1)=%d", mumps->id.INFOG(1));
1196397b6df1SKris Buschelman 
1197b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
11981baa6e33SBarry Smith   if (second_solve) PetscCall(MatMumpsHandleSchur_Private(A, PETSC_TRUE));
1199b5fa320bSStefano Zampini 
12002d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) { /* convert mumps distributed solution to petsc mpi x */
1201a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
1202a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
12039566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&mumps->scat_sol));
1204397b6df1SKris Buschelman     }
1205a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
1206a6053eceSJunchao Zhang       PetscInt *isol2_loc = NULL;
12079566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, mumps->id.lsol_loc, 0, 1, &is_iden)); /* from */
12089566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mumps->id.lsol_loc, &isol2_loc));
1209a6053eceSJunchao Zhang       for (i = 0; i < mumps->id.lsol_loc; i++) isol2_loc[i] = mumps->id.isol_loc[i] - 1;                        /* change Fortran style to C style */
12109566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, mumps->id.lsol_loc, isol2_loc, PETSC_OWN_POINTER, &is_petsc)); /* to */
12119566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(mumps->x_seq, is_iden, x, is_petsc, &mumps->scat_sol));
12129566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_iden));
12139566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_petsc));
1214a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
1215397b6df1SKris Buschelman     }
1216a5e57a09SHong Zhang 
12179566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mumps->scat_sol, mumps->x_seq, x, INSERT_VALUES, SCATTER_FORWARD));
12189566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mumps->scat_sol, mumps->x_seq, x, INSERT_VALUES, SCATTER_FORWARD));
1219329ec9b3SHong Zhang   }
1220353d7d71SJunchao Zhang 
122167602552SJunchao Zhang   if (mumps->petsc_size > 1) {
122225aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
12239566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(b, &rarray));
122425aac85cSJunchao Zhang     } else if (!mumps->myid) {
12259566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mumps->b_seq, &array));
122625aac85cSJunchao Zhang     }
12279566063dSJacob Faibussowitsch   } else PetscCall(VecRestoreArray(x, &array));
1228353d7d71SJunchao Zhang 
12299566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * mumps->id.RINFO(3)));
12303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1231397b6df1SKris Buschelman }
1232397b6df1SKris Buschelman 
1233d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_MUMPS(Mat A, Vec b, Vec x)
1234d71ae5a4SJacob Faibussowitsch {
1235e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)A->data;
123651d5961aSHong Zhang 
123751d5961aSHong Zhang   PetscFunctionBegin;
1238a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
12399566063dSJacob Faibussowitsch   PetscCall(MatSolve_MUMPS(A, b, x));
1240a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
12413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124251d5961aSHong Zhang }
124351d5961aSHong Zhang 
1244d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatSolve_MUMPS(Mat A, Mat B, Mat X)
1245d71ae5a4SJacob Faibussowitsch {
1246b8491c3eSStefano Zampini   Mat                Bt = NULL;
1247a6053eceSJunchao Zhang   PetscBool          denseX, denseB, flg, flgT;
1248e69c285eSBarry Smith   Mat_MUMPS         *mumps = (Mat_MUMPS *)A->data;
1249334c5f61SHong Zhang   PetscInt           i, nrhs, M;
12501683a169SBarry Smith   PetscScalar       *array;
12511683a169SBarry Smith   const PetscScalar *rbray;
1252a6053eceSJunchao Zhang   PetscInt           lsol_loc, nlsol_loc, *idxx, iidx = 0;
1253a6053eceSJunchao Zhang   PetscMUMPSInt     *isol_loc, *isol_loc_save;
12541683a169SBarry Smith   PetscScalar       *bray, *sol_loc, *sol_loc_save;
1255be818407SHong Zhang   IS                 is_to, is_from;
1256beae5ec0SHong Zhang   PetscInt           k, proc, j, m, myrstart;
1257be818407SHong Zhang   const PetscInt    *rstart;
125867602552SJunchao Zhang   Vec                v_mpi, msol_loc;
125967602552SJunchao Zhang   VecScatter         scat_sol;
126067602552SJunchao Zhang   Vec                b_seq;
126167602552SJunchao Zhang   VecScatter         scat_rhs;
1262be818407SHong Zhang   PetscScalar       *aa;
1263be818407SHong Zhang   PetscInt           spnr, *ia, *ja;
1264d56c302dSHong Zhang   Mat_MPIAIJ        *b = NULL;
1265bda8bf91SBarry Smith 
1266e0b74bf9SHong Zhang   PetscFunctionBegin;
12679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &denseX, MATSEQDENSE, MATMPIDENSE, NULL));
126828b400f6SJacob Faibussowitsch   PetscCheck(denseX, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix");
1269be818407SHong Zhang 
12709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &denseB, MATSEQDENSE, MATMPIDENSE, NULL));
1271a6053eceSJunchao Zhang   if (denseB) {
127208401ef6SPierre Jolivet     PetscCheck(B->rmap->n == X->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix B and X must have same row distribution");
1273be818407SHong Zhang     mumps->id.ICNTL(20) = 0; /* dense RHS */
12740e6b8875SHong Zhang   } else {                   /* sparse B */
127508401ef6SPierre Jolivet     PetscCheck(X != B, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_IDN, "X and B must be different matrices");
1276013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &flgT));
1277da81f932SPierre Jolivet     if (flgT) { /* input B is transpose of actual RHS matrix,
12780e6b8875SHong Zhang                  because mumps requires sparse compressed COLUMN storage! See MatMatTransposeSolve_MUMPS() */
12799566063dSJacob Faibussowitsch       PetscCall(MatTransposeGetMat(B, &Bt));
1280013e2dc7SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONG, "Matrix B must be MATTRANSPOSEVIRTUAL matrix");
1281be818407SHong Zhang     mumps->id.ICNTL(20) = 1; /* sparse RHS */
1282b8491c3eSStefano Zampini   }
128387b22cf4SHong Zhang 
12849566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, &nrhs));
12859481e6e9SHong Zhang   mumps->id.nrhs = nrhs;
12869481e6e9SHong Zhang   mumps->id.lrhs = M;
12872b691707SHong Zhang   mumps->id.rhs  = NULL;
12889481e6e9SHong Zhang 
12892d4298aeSJunchao Zhang   if (mumps->petsc_size == 1) {
1290b8491c3eSStefano Zampini     PetscScalar *aa;
1291b8491c3eSStefano Zampini     PetscInt     spnr, *ia, *ja;
1292e94cce23SStefano Zampini     PetscBool    second_solve = PETSC_FALSE;
1293b8491c3eSStefano Zampini 
12949566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(X, &array));
1295b8491c3eSStefano Zampini     mumps->id.rhs = (MumpsScalar *)array;
12962b691707SHong Zhang 
1297a6053eceSJunchao Zhang     if (denseB) {
12982b691707SHong Zhang       /* copy B to X */
12999566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(B, &rbray));
13009566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array, rbray, M * nrhs));
13019566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(B, &rbray));
13022b691707SHong Zhang     } else { /* sparse B */
13039566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(Bt, &aa));
13049566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Bt, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
130528b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
13069566063dSJacob Faibussowitsch       PetscCall(PetscMUMPSIntCSRCast(mumps, spnr, ia, ja, &mumps->id.irhs_ptr, &mumps->id.irhs_sparse, &mumps->id.nz_rhs));
1307b8491c3eSStefano Zampini       mumps->id.rhs_sparse = (MumpsScalar *)aa;
1308b8491c3eSStefano Zampini     }
1309e94cce23SStefano Zampini     /* handle condensation step of Schur complement (if any) */
1310583f777eSStefano Zampini     if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
1311e94cce23SStefano Zampini       second_solve = PETSC_TRUE;
13129566063dSJacob Faibussowitsch       PetscCall(MatMumpsHandleSchur_Private(A, PETSC_FALSE));
1313e94cce23SStefano Zampini     }
13142cd7d884SHong Zhang     /* solve phase */
13152cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
13163ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
131708401ef6SPierre Jolivet     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in solve phase: INFOG(1)=%d", mumps->id.INFOG(1));
1318b5fa320bSStefano Zampini 
1319b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
13201baa6e33SBarry Smith     if (second_solve) PetscCall(MatMumpsHandleSchur_Private(A, PETSC_TRUE));
1321a6053eceSJunchao Zhang     if (!denseB) { /* sparse B */
13229566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(Bt, &aa));
13239566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Bt, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
132428b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot restore IJ structure");
1325b8491c3eSStefano Zampini     }
13269566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X, &array));
13273ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1328be818407SHong Zhang   }
1329801fbe65SHong Zhang 
13302ef1f0ffSBarry Smith   /* parallel case: MUMPS requires rhs B to be centralized on the host! */
1331aed4548fSBarry Smith   PetscCheck(mumps->petsc_size <= 1 || !mumps->id.ICNTL(19), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Parallel Schur complements not yet supported from PETSc");
1332241dbb5eSStefano Zampini 
1333beae5ec0SHong Zhang   /* create msol_loc to hold mumps local solution */
13341683a169SBarry Smith   isol_loc_save = mumps->id.isol_loc; /* save it for MatSolve() */
13351683a169SBarry Smith   sol_loc_save  = (PetscScalar *)mumps->id.sol_loc;
1336801fbe65SHong Zhang 
1337a1dfcbd9SJunchao Zhang   lsol_loc  = mumps->id.lsol_loc;
133871aed81dSHong Zhang   nlsol_loc = nrhs * lsol_loc; /* length of sol_loc */
13399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nlsol_loc, &sol_loc, lsol_loc, &isol_loc));
1340940cd9d6SSatish Balay   mumps->id.sol_loc  = (MumpsScalar *)sol_loc;
1341801fbe65SHong Zhang   mumps->id.isol_loc = isol_loc;
1342801fbe65SHong Zhang 
13439566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nlsol_loc, (PetscScalar *)sol_loc, &msol_loc));
13442cd7d884SHong Zhang 
134567602552SJunchao Zhang   if (denseB) {
134625aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
134767602552SJunchao Zhang       mumps->id.ICNTL(20) = 10; /* dense distributed RHS */
13489566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(B, &rbray));
13499566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetUpDistRHSInfo(A, nrhs, rbray));
13509566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(B, &rbray));
13519566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, &m, NULL));
13529566063dSJacob Faibussowitsch       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)B), 1, nrhs * m, nrhs * M, NULL, &v_mpi));
135325aac85cSJunchao Zhang     } else {
135425aac85cSJunchao Zhang       mumps->id.ICNTL(20) = 0; /* dense centralized RHS */
135580577c12SJunchao Zhang       /* TODO: Because of non-contiguous indices, the created vecscatter scat_rhs is not done in MPI_Gather, resulting in
135680577c12SJunchao Zhang         very inefficient communication. An optimization is to use VecScatterCreateToZero to gather B to rank 0. Then on rank
135780577c12SJunchao Zhang         0, re-arrange B into desired order, which is a local operation.
135880577c12SJunchao Zhang       */
135980577c12SJunchao Zhang 
136067602552SJunchao Zhang       /* scatter v_mpi to b_seq because MUMPS before 5.3.0 only supports centralized rhs */
1361be818407SHong Zhang       /* wrap dense rhs matrix B into a vector v_mpi */
13629566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, &m, NULL));
13639566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArray(B, &bray));
13649566063dSJacob Faibussowitsch       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)B), 1, nrhs * m, nrhs * M, (const PetscScalar *)bray, &v_mpi));
13659566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArray(B, &bray));
13662b691707SHong Zhang 
1367be818407SHong Zhang       /* scatter v_mpi to b_seq in proc[0]. MUMPS requires rhs to be centralized on the host! */
1368801fbe65SHong Zhang       if (!mumps->myid) {
1369beae5ec0SHong Zhang         PetscInt *idx;
1370beae5ec0SHong Zhang         /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B */
13719566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nrhs * M, &idx));
13729566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRanges(B, &rstart));
1373be818407SHong Zhang         k = 0;
13742d4298aeSJunchao Zhang         for (proc = 0; proc < mumps->petsc_size; proc++) {
1375be818407SHong Zhang           for (j = 0; j < nrhs; j++) {
1376beae5ec0SHong Zhang             for (i = rstart[proc]; i < rstart[proc + 1]; i++) idx[k++] = j * M + i;
1377be818407SHong Zhang           }
1378be818407SHong Zhang         }
1379be818407SHong Zhang 
13809566063dSJacob Faibussowitsch         PetscCall(VecCreateSeq(PETSC_COMM_SELF, nrhs * M, &b_seq));
13819566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrhs * M, idx, PETSC_OWN_POINTER, &is_to));
13829566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF, nrhs * M, 0, 1, &is_from));
1383801fbe65SHong Zhang       } else {
13849566063dSJacob Faibussowitsch         PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &b_seq));
13859566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is_to));
13869566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is_from));
1387801fbe65SHong Zhang       }
13889566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(v_mpi, is_from, b_seq, is_to, &scat_rhs));
13899566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(scat_rhs, v_mpi, b_seq, INSERT_VALUES, SCATTER_FORWARD));
13909566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_to));
13919566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_from));
13929566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(scat_rhs, v_mpi, b_seq, INSERT_VALUES, SCATTER_FORWARD));
1393801fbe65SHong Zhang 
1394801fbe65SHong Zhang       if (!mumps->myid) { /* define rhs on the host */
13959566063dSJacob Faibussowitsch         PetscCall(VecGetArray(b_seq, &bray));
1396940cd9d6SSatish Balay         mumps->id.rhs = (MumpsScalar *)bray;
13979566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(b_seq, &bray));
1398801fbe65SHong Zhang       }
139925aac85cSJunchao Zhang     }
14002b691707SHong Zhang   } else { /* sparse B */
14012b691707SHong Zhang     b = (Mat_MPIAIJ *)Bt->data;
14022b691707SHong Zhang 
1403be818407SHong Zhang     /* wrap dense X into a vector v_mpi */
14049566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(X, &m, NULL));
14059566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(X, &bray));
14069566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)X), 1, nrhs * m, nrhs * M, (const PetscScalar *)bray, &v_mpi));
14079566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X, &bray));
14082b691707SHong Zhang 
14092b691707SHong Zhang     if (!mumps->myid) {
14109566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(b->A, &aa));
14119566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(b->A, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
141228b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
14139566063dSJacob Faibussowitsch       PetscCall(PetscMUMPSIntCSRCast(mumps, spnr, ia, ja, &mumps->id.irhs_ptr, &mumps->id.irhs_sparse, &mumps->id.nz_rhs));
14142b691707SHong Zhang       mumps->id.rhs_sparse = (MumpsScalar *)aa;
14152b691707SHong Zhang     } else {
14162b691707SHong Zhang       mumps->id.irhs_ptr    = NULL;
14172b691707SHong Zhang       mumps->id.irhs_sparse = NULL;
14182b691707SHong Zhang       mumps->id.nz_rhs      = 0;
14192b691707SHong Zhang       mumps->id.rhs_sparse  = NULL;
14202b691707SHong Zhang     }
14212b691707SHong Zhang   }
14222b691707SHong Zhang 
1423801fbe65SHong Zhang   /* solve phase */
1424801fbe65SHong Zhang   mumps->id.job = JOB_SOLVE;
14253ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
142608401ef6SPierre Jolivet   PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in solve phase: INFOG(1)=%d", mumps->id.INFOG(1));
1427801fbe65SHong Zhang 
1428334c5f61SHong Zhang   /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
14299566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(X, &array));
14309566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(v_mpi, array));
1431801fbe65SHong Zhang 
1432334c5f61SHong Zhang   /* create scatter scat_sol */
14339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(X, &rstart));
1434beae5ec0SHong Zhang   /* iidx: index for scatter mumps solution to petsc X */
1435beae5ec0SHong Zhang 
14369566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, nlsol_loc, 0, 1, &is_from));
14379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nlsol_loc, &idxx));
1438beae5ec0SHong Zhang   for (i = 0; i < lsol_loc; i++) {
1439beae5ec0SHong Zhang     isol_loc[i] -= 1; /* change Fortran style to C style. isol_loc[i+j*lsol_loc] contains x[isol_loc[i]] in j-th vector */
1440beae5ec0SHong Zhang 
14412d4298aeSJunchao Zhang     for (proc = 0; proc < mumps->petsc_size; proc++) {
1442beae5ec0SHong Zhang       if (isol_loc[i] >= rstart[proc] && isol_loc[i] < rstart[proc + 1]) {
1443beae5ec0SHong Zhang         myrstart = rstart[proc];
1444beae5ec0SHong Zhang         k        = isol_loc[i] - myrstart;          /* local index on 1st column of petsc vector X */
1445beae5ec0SHong Zhang         iidx     = k + myrstart * nrhs;             /* maps mumps isol_loc[i] to petsc index in X */
1446beae5ec0SHong Zhang         m        = rstart[proc + 1] - rstart[proc]; /* rows of X for this proc */
1447beae5ec0SHong Zhang         break;
1448be818407SHong Zhang       }
1449be818407SHong Zhang     }
1450be818407SHong Zhang 
1451beae5ec0SHong Zhang     for (j = 0; j < nrhs; j++) idxx[i + j * lsol_loc] = iidx + j * m;
1452801fbe65SHong Zhang   }
14539566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nlsol_loc, idxx, PETSC_COPY_VALUES, &is_to));
14549566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(msol_loc, is_from, v_mpi, is_to, &scat_sol));
14559566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat_sol, msol_loc, v_mpi, INSERT_VALUES, SCATTER_FORWARD));
14569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_from));
14579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_to));
14589566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat_sol, msol_loc, v_mpi, INSERT_VALUES, SCATTER_FORWARD));
14599566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(X, &array));
146071aed81dSHong Zhang 
146171aed81dSHong Zhang   /* free spaces */
14621683a169SBarry Smith   mumps->id.sol_loc  = (MumpsScalar *)sol_loc_save;
146371aed81dSHong Zhang   mumps->id.isol_loc = isol_loc_save;
146471aed81dSHong Zhang 
14659566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sol_loc, isol_loc));
14669566063dSJacob Faibussowitsch   PetscCall(PetscFree(idxx));
14679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&msol_loc));
14689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v_mpi));
1469a6053eceSJunchao Zhang   if (!denseB) {
14702b691707SHong Zhang     if (!mumps->myid) {
1471d56c302dSHong Zhang       b = (Mat_MPIAIJ *)Bt->data;
14729566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(b->A, &aa));
14739566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(b->A, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
147428b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot restore IJ structure");
14752b691707SHong Zhang     }
14762b691707SHong Zhang   } else {
147725aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) {
14789566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b_seq));
14799566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&scat_rhs));
148025aac85cSJunchao Zhang     }
14812b691707SHong Zhang   }
14829566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scat_sol));
14839566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * nrhs * mumps->id.RINFO(3)));
14843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1485e0b74bf9SHong Zhang }
1486e0b74bf9SHong Zhang 
1487d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatSolveTranspose_MUMPS(Mat A, Mat B, Mat X)
1488d71ae5a4SJacob Faibussowitsch {
1489b18964edSHong Zhang   Mat_MUMPS    *mumps    = (Mat_MUMPS *)A->data;
1490b18964edSHong Zhang   PetscMUMPSInt oldvalue = mumps->id.ICNTL(9);
1491b18964edSHong Zhang 
1492b18964edSHong Zhang   PetscFunctionBegin;
1493b18964edSHong Zhang   mumps->id.ICNTL(9) = 0;
1494b18964edSHong Zhang   PetscCall(MatMatSolve_MUMPS(A, B, X));
1495b18964edSHong Zhang   mumps->id.ICNTL(9) = oldvalue;
14963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1497b18964edSHong Zhang }
1498b18964edSHong Zhang 
1499d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatTransposeSolve_MUMPS(Mat A, Mat Bt, Mat X)
1500d71ae5a4SJacob Faibussowitsch {
1501eb3ef3b2SHong Zhang   PetscBool flg;
1502eb3ef3b2SHong Zhang   Mat       B;
1503eb3ef3b2SHong Zhang 
1504eb3ef3b2SHong Zhang   PetscFunctionBegin;
15059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)Bt, &flg, MATSEQAIJ, MATMPIAIJ, NULL));
150628b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)Bt), PETSC_ERR_ARG_WRONG, "Matrix Bt must be MATAIJ matrix");
1507eb3ef3b2SHong Zhang 
1508eb3ef3b2SHong Zhang   /* Create B=Bt^T that uses Bt's data structure */
15099566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(Bt, &B));
1510eb3ef3b2SHong Zhang 
15119566063dSJacob Faibussowitsch   PetscCall(MatMatSolve_MUMPS(A, B, X));
15129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
15133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1514eb3ef3b2SHong Zhang }
1515eb3ef3b2SHong Zhang 
1516ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1517a58c3f20SHong Zhang /*
1518a58c3f20SHong Zhang   input:
1519a58c3f20SHong Zhang    F:        numeric factor
1520a58c3f20SHong Zhang   output:
1521a58c3f20SHong Zhang    nneg:     total number of negative pivots
152219d49a3bSHong Zhang    nzero:    total number of zero pivots
152319d49a3bSHong Zhang    npos:     (global dimension of F) - nneg - nzero
1524a58c3f20SHong Zhang */
1525d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F, PetscInt *nneg, PetscInt *nzero, PetscInt *npos)
1526d71ae5a4SJacob Faibussowitsch {
1527e69c285eSBarry Smith   Mat_MUMPS  *mumps = (Mat_MUMPS *)F->data;
1528c1490034SHong Zhang   PetscMPIInt size;
1529a58c3f20SHong Zhang 
1530a58c3f20SHong Zhang   PetscFunctionBegin;
15319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)F), &size));
1532bcb30aebSHong Zhang   /* MUMPS 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK */
1533aed4548fSBarry Smith   PetscCheck(size <= 1 || mumps->id.ICNTL(13) == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia", mumps->id.INFOG(13));
1534ed85ac9fSHong Zhang 
1535710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1536ed85ac9fSHong Zhang   if (nzero || npos) {
153708401ef6SPierre Jolivet     PetscCheck(mumps->id.ICNTL(24) == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
1538710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1539710ac8efSHong Zhang     if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1540a58c3f20SHong Zhang   }
15413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1542a58c3f20SHong Zhang }
154319d49a3bSHong Zhang #endif
1544a58c3f20SHong Zhang 
1545d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGatherNonzerosOnMaster(MatReuse reuse, Mat_MUMPS *mumps)
1546d71ae5a4SJacob Faibussowitsch {
1547a6053eceSJunchao Zhang   PetscInt       i, nreqs;
1548a6053eceSJunchao Zhang   PetscMUMPSInt *irn, *jcn;
1549a6053eceSJunchao Zhang   PetscMPIInt    count;
1550a6053eceSJunchao Zhang   PetscInt64     totnnz, remain;
1551a6053eceSJunchao Zhang   const PetscInt osize = mumps->omp_comm_size;
1552a6053eceSJunchao Zhang   PetscScalar   *val;
15533ab56b82SJunchao Zhang 
15543ab56b82SJunchao Zhang   PetscFunctionBegin;
1555a6053eceSJunchao Zhang   if (osize > 1) {
15563ab56b82SJunchao Zhang     if (reuse == MAT_INITIAL_MATRIX) {
15573ab56b82SJunchao Zhang       /* master first gathers counts of nonzeros to receive */
15589566063dSJacob Faibussowitsch       if (mumps->is_omp_master) PetscCall(PetscMalloc1(osize, &mumps->recvcount));
15599566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gather(&mumps->nnz, 1, MPIU_INT64, mumps->recvcount, 1, MPIU_INT64, 0 /*master*/, mumps->omp_comm));
15603ab56b82SJunchao Zhang 
1561a6053eceSJunchao Zhang       /* Then each computes number of send/recvs */
15623ab56b82SJunchao Zhang       if (mumps->is_omp_master) {
1563a6053eceSJunchao Zhang         /* Start from 1 since self communication is not done in MPI */
1564a6053eceSJunchao Zhang         nreqs = 0;
1565a6053eceSJunchao Zhang         for (i = 1; i < osize; i++) nreqs += (mumps->recvcount[i] + PETSC_MPI_INT_MAX - 1) / PETSC_MPI_INT_MAX;
1566a6053eceSJunchao Zhang       } else {
1567a6053eceSJunchao Zhang         nreqs = (mumps->nnz + PETSC_MPI_INT_MAX - 1) / PETSC_MPI_INT_MAX;
15683ab56b82SJunchao Zhang       }
156935cb6cd3SPierre Jolivet       PetscCall(PetscMalloc1(nreqs * 3, &mumps->reqs)); /* Triple the requests since we send irn, jcn and val separately */
15703ab56b82SJunchao Zhang 
1571a6053eceSJunchao Zhang       /* The following code is doing a very simple thing: omp_master rank gathers irn/jcn/val from others.
1572a6053eceSJunchao Zhang          MPI_Gatherv would be enough if it supports big counts > 2^31-1. Since it does not, and mumps->nnz
1573a6053eceSJunchao Zhang          might be a prime number > 2^31-1, we have to slice the message. Note omp_comm_size
1574a6053eceSJunchao Zhang          is very small, the current approach should have no extra overhead compared to MPI_Gatherv.
1575a6053eceSJunchao Zhang        */
1576a6053eceSJunchao Zhang       nreqs = 0; /* counter for actual send/recvs */
15773ab56b82SJunchao Zhang       if (mumps->is_omp_master) {
1578a6053eceSJunchao Zhang         for (i = 0, totnnz = 0; i < osize; i++) totnnz += mumps->recvcount[i]; /* totnnz = sum of nnz over omp_comm */
15799566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(totnnz, &irn, totnnz, &jcn));
15809566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(totnnz, &val));
1581a6053eceSJunchao Zhang 
1582a6053eceSJunchao Zhang         /* Self communication */
15839566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(irn, mumps->irn, mumps->nnz));
15849566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jcn, mumps->jcn, mumps->nnz));
15859566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(val, mumps->val, mumps->nnz));
1586a6053eceSJunchao Zhang 
1587a6053eceSJunchao Zhang         /* Replace mumps->irn/jcn etc on master with the newly allocated bigger arrays */
15889566063dSJacob Faibussowitsch         PetscCall(PetscFree2(mumps->irn, mumps->jcn));
15899566063dSJacob Faibussowitsch         PetscCall(PetscFree(mumps->val_alloc));
1590a6053eceSJunchao Zhang         mumps->nnz = totnnz;
15913ab56b82SJunchao Zhang         mumps->irn = irn;
15923ab56b82SJunchao Zhang         mumps->jcn = jcn;
1593a6053eceSJunchao Zhang         mumps->val = mumps->val_alloc = val;
1594a6053eceSJunchao Zhang 
1595a6053eceSJunchao Zhang         irn += mumps->recvcount[0]; /* recvcount[0] is old mumps->nnz on omp rank 0 */
1596a6053eceSJunchao Zhang         jcn += mumps->recvcount[0];
1597a6053eceSJunchao Zhang         val += mumps->recvcount[0];
1598a6053eceSJunchao Zhang 
1599a6053eceSJunchao Zhang         /* Remote communication */
1600a6053eceSJunchao Zhang         for (i = 1; i < osize; i++) {
1601a6053eceSJunchao Zhang           count  = PetscMin(mumps->recvcount[i], PETSC_MPI_INT_MAX);
1602a6053eceSJunchao Zhang           remain = mumps->recvcount[i] - count;
1603a6053eceSJunchao Zhang           while (count > 0) {
16049566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(irn, count, MPIU_MUMPSINT, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
16059566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(jcn, count, MPIU_MUMPSINT, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
16069566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(val, count, MPIU_SCALAR, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1607a6053eceSJunchao Zhang             irn += count;
1608a6053eceSJunchao Zhang             jcn += count;
1609a6053eceSJunchao Zhang             val += count;
1610a6053eceSJunchao Zhang             count = PetscMin(remain, PETSC_MPI_INT_MAX);
1611a6053eceSJunchao Zhang             remain -= count;
1612a6053eceSJunchao Zhang           }
16133ab56b82SJunchao Zhang         }
16143ab56b82SJunchao Zhang       } else {
1615a6053eceSJunchao Zhang         irn    = mumps->irn;
1616a6053eceSJunchao Zhang         jcn    = mumps->jcn;
1617a6053eceSJunchao Zhang         val    = mumps->val;
1618a6053eceSJunchao Zhang         count  = PetscMin(mumps->nnz, PETSC_MPI_INT_MAX);
1619a6053eceSJunchao Zhang         remain = mumps->nnz - count;
1620a6053eceSJunchao Zhang         while (count > 0) {
16219566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(irn, count, MPIU_MUMPSINT, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
16229566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(jcn, count, MPIU_MUMPSINT, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
16239566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(val, count, MPIU_SCALAR, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1624a6053eceSJunchao Zhang           irn += count;
1625a6053eceSJunchao Zhang           jcn += count;
1626a6053eceSJunchao Zhang           val += count;
1627a6053eceSJunchao Zhang           count = PetscMin(remain, PETSC_MPI_INT_MAX);
1628a6053eceSJunchao Zhang           remain -= count;
16293ab56b82SJunchao Zhang         }
16303ab56b82SJunchao Zhang       }
1631a6053eceSJunchao Zhang     } else {
1632a6053eceSJunchao Zhang       nreqs = 0;
1633a6053eceSJunchao Zhang       if (mumps->is_omp_master) {
1634a6053eceSJunchao Zhang         val = mumps->val + mumps->recvcount[0];
1635a6053eceSJunchao Zhang         for (i = 1; i < osize; i++) { /* Remote communication only since self data is already in place */
1636a6053eceSJunchao Zhang           count  = PetscMin(mumps->recvcount[i], PETSC_MPI_INT_MAX);
1637a6053eceSJunchao Zhang           remain = mumps->recvcount[i] - count;
1638a6053eceSJunchao Zhang           while (count > 0) {
16399566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(val, count, MPIU_SCALAR, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1640a6053eceSJunchao Zhang             val += count;
1641a6053eceSJunchao Zhang             count = PetscMin(remain, PETSC_MPI_INT_MAX);
1642a6053eceSJunchao Zhang             remain -= count;
1643a6053eceSJunchao Zhang           }
1644a6053eceSJunchao Zhang         }
1645a6053eceSJunchao Zhang       } else {
1646a6053eceSJunchao Zhang         val    = mumps->val;
1647a6053eceSJunchao Zhang         count  = PetscMin(mumps->nnz, PETSC_MPI_INT_MAX);
1648a6053eceSJunchao Zhang         remain = mumps->nnz - count;
1649a6053eceSJunchao Zhang         while (count > 0) {
16509566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(val, count, MPIU_SCALAR, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1651a6053eceSJunchao Zhang           val += count;
1652a6053eceSJunchao Zhang           count = PetscMin(remain, PETSC_MPI_INT_MAX);
1653a6053eceSJunchao Zhang           remain -= count;
1654a6053eceSJunchao Zhang         }
1655a6053eceSJunchao Zhang       }
1656a6053eceSJunchao Zhang     }
16579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nreqs, mumps->reqs, MPI_STATUSES_IGNORE));
1658a6053eceSJunchao Zhang     mumps->tag++; /* It is totally fine for above send/recvs to share one mpi tag */
1659a6053eceSJunchao Zhang   }
16603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16613ab56b82SJunchao Zhang }
16623ab56b82SJunchao Zhang 
1663d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorNumeric_MUMPS(Mat F, Mat A, const MatFactorInfo *info)
1664d71ae5a4SJacob Faibussowitsch {
1665e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)(F)->data;
1666ace3abfcSBarry Smith   PetscBool  isMPIAIJ;
1667397b6df1SKris Buschelman 
1668397b6df1SKris Buschelman   PetscFunctionBegin;
1669dbf6bb8dSprj-   if (mumps->id.INFOG(1) < 0 && !(mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0)) {
167048a46eb9SPierre Jolivet     if (mumps->id.INFOG(1) == -6) PetscCall(PetscInfo(A, "MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
16719566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
16723ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
16732aca8efcSHong Zhang   }
16746baea169SHong Zhang 
16759566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, mumps));
16769566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_REUSE_MATRIX, mumps));
1677397b6df1SKris Buschelman 
1678397b6df1SKris Buschelman   /* numerical factorization phase */
1679a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
16804e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1681ad540459SPierre Jolivet     if (!mumps->myid) mumps->id.a = (MumpsScalar *)mumps->val;
1682397b6df1SKris Buschelman   } else {
1683940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar *)mumps->val;
1684397b6df1SKris Buschelman   }
16853ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
1686a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
16877a46b595SBarry Smith     PetscCheck(!A->erroriffailure, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d", mumps->id.INFOG(1), mumps->id.INFO(2));
1688c0d63f2fSHong Zhang     if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
16899566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1690603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1691c0d63f2fSHong Zhang     } else if (mumps->id.INFOG(1) == -13) {
16929566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1693603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1694c0d63f2fSHong Zhang     } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10)) {
16959566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1696603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OUTMEMORY;
16972aca8efcSHong Zhang     } else {
16989566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1699603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OTHER;
1700151787a6SHong Zhang     }
17012aca8efcSHong Zhang   }
1702aed4548fSBarry Smith   PetscCheck(mumps->myid || mumps->id.ICNTL(16) <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "  mumps->id.ICNTL(16):=%d", mumps->id.INFOG(16));
1703397b6df1SKris Buschelman 
1704b3cb21ddSStefano Zampini   F->assembled = PETSC_TRUE;
1705d47f36abSHong Zhang 
1706b3cb21ddSStefano Zampini   if (F->schur) { /* reset Schur status to unfactored */
17073cb7dd0eSStefano Zampini #if defined(PETSC_HAVE_CUDA)
1708c70f7ee4SJunchao Zhang     F->schur->offloadmask = PETSC_OFFLOAD_CPU;
17093cb7dd0eSStefano Zampini #endif
1710b3cb21ddSStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1711b3cb21ddSStefano Zampini       mumps->id.ICNTL(19) = 2;
17129566063dSJacob Faibussowitsch       PetscCall(MatTranspose(F->schur, MAT_INPLACE_MATRIX, &F->schur));
1713b3cb21ddSStefano Zampini     }
17149566063dSJacob Faibussowitsch     PetscCall(MatFactorRestoreSchurComplement(F, NULL, MAT_FACTOR_SCHUR_UNFACTORED));
1715b3cb21ddSStefano Zampini   }
171667877ebaSShri Abhyankar 
1717066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1718066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1719066565c5SStefano Zampini 
17203ab56b82SJunchao Zhang   if (!mumps->is_omp_master) mumps->id.INFO(23) = 0;
17212d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
172267877ebaSShri Abhyankar     PetscInt     lsol_loc;
172367877ebaSShri Abhyankar     PetscScalar *sol_loc;
17242205254eSKarl Rupp 
17259566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ));
1726c2093ab7SHong Zhang 
1727c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
1728c2093ab7SHong Zhang     if (mumps->x_seq) {
17299566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&mumps->scat_sol));
17309566063dSJacob Faibussowitsch       PetscCall(PetscFree2(mumps->id.sol_loc, mumps->id.isol_loc));
17319566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mumps->x_seq));
1732c2093ab7SHong Zhang     }
1733a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
17349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lsol_loc, &sol_loc, lsol_loc, &mumps->id.isol_loc));
1735a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
1736940cd9d6SSatish Balay     mumps->id.sol_loc  = (MumpsScalar *)sol_loc;
17379566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, lsol_loc, sol_loc, &mumps->x_seq));
173867877ebaSShri Abhyankar   }
17399566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(mumps->id.RINFO(2)));
17403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1741397b6df1SKris Buschelman }
1742397b6df1SKris Buschelman 
17439a2535b5SHong Zhang /* Sets MUMPS options from the options database */
1744d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetFromOptions_MUMPS(Mat F, Mat A)
1745d71ae5a4SJacob Faibussowitsch {
1746e69c285eSBarry Smith   Mat_MUMPS    *mumps = (Mat_MUMPS *)F->data;
1747413bcc21SPierre Jolivet   PetscMUMPSInt icntl = 0, size, *listvar_schur;
174845e3843bSPierre Jolivet   PetscInt      info[80], i, ninfo = 80, rbs, cbs;
1749413bcc21SPierre Jolivet   PetscBool     flg = PETSC_FALSE, schur = (PetscBool)(mumps->id.ICNTL(26) == -1);
1750413bcc21SPierre Jolivet   MumpsScalar  *arr;
1751dcd589f8SShri Abhyankar 
1752dcd589f8SShri Abhyankar   PetscFunctionBegin;
175326cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "MUMPS Options", "Mat");
1754413bcc21SPierre Jolivet   if (mumps->id.job == JOB_NULL) { /* MatSetFromOptions_MUMPS() has never been called before */
1755413bcc21SPierre Jolivet     PetscInt nthreads   = 0;
1756413bcc21SPierre Jolivet     PetscInt nCNTL_pre  = mumps->CNTL_pre ? mumps->CNTL_pre[0] : 0;
1757413bcc21SPierre Jolivet     PetscInt nICNTL_pre = mumps->ICNTL_pre ? mumps->ICNTL_pre[0] : 0;
1758413bcc21SPierre Jolivet 
1759413bcc21SPierre Jolivet     mumps->petsc_comm = PetscObjectComm((PetscObject)A);
1760413bcc21SPierre Jolivet     PetscCallMPI(MPI_Comm_size(mumps->petsc_comm, &mumps->petsc_size));
1761413bcc21SPierre Jolivet     PetscCallMPI(MPI_Comm_rank(mumps->petsc_comm, &mumps->myid)); /* "if (!myid)" still works even if mumps_comm is different */
1762413bcc21SPierre Jolivet 
1763413bcc21SPierre Jolivet     PetscCall(PetscOptionsName("-mat_mumps_use_omp_threads", "Convert MPI processes into OpenMP threads", "None", &mumps->use_petsc_omp_support));
1764413bcc21SPierre Jolivet     if (mumps->use_petsc_omp_support) nthreads = -1; /* -1 will let PetscOmpCtrlCreate() guess a proper value when user did not supply one */
1765413bcc21SPierre Jolivet     /* do not use PetscOptionsInt() so that the option -mat_mumps_use_omp_threads is not displayed twice in the help */
1766413bcc21SPierre Jolivet     PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)F)->prefix, "-mat_mumps_use_omp_threads", &nthreads, NULL));
1767413bcc21SPierre Jolivet     if (mumps->use_petsc_omp_support) {
17689371c9d4SSatish Balay       PetscCheck(PetscDefined(HAVE_OPENMP_SUPPORT), PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "The system does not have PETSc OpenMP support but you added the -%smat_mumps_use_omp_threads option. Configure PETSc with --with-openmp --download-hwloc (or --with-hwloc) to enable it, see more in MATSOLVERMUMPS manual",
17699371c9d4SSatish Balay                  ((PetscObject)F)->prefix ? ((PetscObject)F)->prefix : "");
1770413bcc21SPierre Jolivet       PetscCheck(!schur, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use -%smat_mumps_use_omp_threads with the Schur complement feature", ((PetscObject)F)->prefix ? ((PetscObject)F)->prefix : "");
1771413bcc21SPierre Jolivet #if defined(PETSC_HAVE_OPENMP_SUPPORT)
1772413bcc21SPierre Jolivet       PetscCall(PetscOmpCtrlCreate(mumps->petsc_comm, nthreads, &mumps->omp_ctrl));
1773413bcc21SPierre Jolivet       PetscCall(PetscOmpCtrlGetOmpComms(mumps->omp_ctrl, &mumps->omp_comm, &mumps->mumps_comm, &mumps->is_omp_master));
1774413bcc21SPierre Jolivet #endif
1775413bcc21SPierre Jolivet     } else {
1776413bcc21SPierre Jolivet       mumps->omp_comm      = PETSC_COMM_SELF;
1777413bcc21SPierre Jolivet       mumps->mumps_comm    = mumps->petsc_comm;
1778413bcc21SPierre Jolivet       mumps->is_omp_master = PETSC_TRUE;
1779413bcc21SPierre Jolivet     }
1780413bcc21SPierre Jolivet     PetscCallMPI(MPI_Comm_size(mumps->omp_comm, &mumps->omp_comm_size));
1781413bcc21SPierre Jolivet     mumps->reqs = NULL;
1782413bcc21SPierre Jolivet     mumps->tag  = 0;
1783413bcc21SPierre Jolivet 
1784413bcc21SPierre Jolivet     if (mumps->mumps_comm != MPI_COMM_NULL) {
1785413bcc21SPierre Jolivet       if (PetscDefined(HAVE_OPENMP_SUPPORT) && mumps->use_petsc_omp_support) {
1786413bcc21SPierre Jolivet         /* It looks like MUMPS does not dup the input comm. Dup a new comm for MUMPS to avoid any tag mismatches. */
1787413bcc21SPierre Jolivet         MPI_Comm comm;
1788413bcc21SPierre Jolivet         PetscCallMPI(MPI_Comm_dup(mumps->mumps_comm, &comm));
1789413bcc21SPierre Jolivet         mumps->mumps_comm = comm;
1790413bcc21SPierre Jolivet       } else PetscCall(PetscCommGetComm(mumps->petsc_comm, &mumps->mumps_comm));
1791413bcc21SPierre Jolivet     }
1792413bcc21SPierre Jolivet 
1793413bcc21SPierre Jolivet     mumps->id.comm_fortran = MPI_Comm_c2f(mumps->mumps_comm);
1794413bcc21SPierre Jolivet     mumps->id.job          = JOB_INIT;
1795413bcc21SPierre Jolivet     mumps->id.par          = 1; /* host participates factorizaton and solve */
1796413bcc21SPierre Jolivet     mumps->id.sym          = mumps->sym;
1797413bcc21SPierre Jolivet 
1798413bcc21SPierre Jolivet     size          = mumps->id.size_schur;
1799413bcc21SPierre Jolivet     arr           = mumps->id.schur;
1800413bcc21SPierre Jolivet     listvar_schur = mumps->id.listvar_schur;
1801413bcc21SPierre Jolivet     PetscMUMPS_c(mumps);
1802413bcc21SPierre Jolivet     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS: INFOG(1)=%d", mumps->id.INFOG(1));
1803413bcc21SPierre Jolivet     /* restore cached ICNTL and CNTL values */
1804413bcc21SPierre Jolivet     for (icntl = 0; icntl < nICNTL_pre; ++icntl) mumps->id.ICNTL(mumps->ICNTL_pre[1 + 2 * icntl]) = mumps->ICNTL_pre[2 + 2 * icntl];
1805413bcc21SPierre Jolivet     for (icntl = 0; icntl < nCNTL_pre; ++icntl) mumps->id.CNTL((PetscInt)mumps->CNTL_pre[1 + 2 * icntl]) = mumps->CNTL_pre[2 + 2 * icntl];
1806413bcc21SPierre Jolivet     PetscCall(PetscFree(mumps->ICNTL_pre));
1807413bcc21SPierre Jolivet     PetscCall(PetscFree(mumps->CNTL_pre));
1808413bcc21SPierre Jolivet 
1809413bcc21SPierre Jolivet     if (schur) {
1810413bcc21SPierre Jolivet       mumps->id.size_schur    = size;
1811413bcc21SPierre Jolivet       mumps->id.schur_lld     = size;
1812413bcc21SPierre Jolivet       mumps->id.schur         = arr;
1813413bcc21SPierre Jolivet       mumps->id.listvar_schur = listvar_schur;
1814413bcc21SPierre Jolivet       if (mumps->petsc_size > 1) {
1815413bcc21SPierre Jolivet         PetscBool gs; /* gs is false if any rank other than root has non-empty IS */
1816413bcc21SPierre Jolivet 
1817413bcc21SPierre Jolivet         mumps->id.ICNTL(19) = 1;                                                                            /* MUMPS returns Schur centralized on the host */
1818413bcc21SPierre Jolivet         gs                  = mumps->myid ? (mumps->id.size_schur ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE; /* always true on root; false on others if their size != 0 */
1819712fec58SPierre Jolivet         PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &gs, 1, MPIU_BOOL, MPI_LAND, mumps->petsc_comm));
1820413bcc21SPierre Jolivet         PetscCheck(gs, PETSC_COMM_SELF, PETSC_ERR_SUP, "MUMPS distributed parallel Schur complements not yet supported from PETSc");
1821413bcc21SPierre Jolivet       } else {
1822413bcc21SPierre Jolivet         if (F->factortype == MAT_FACTOR_LU) {
1823413bcc21SPierre Jolivet           mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1824413bcc21SPierre Jolivet         } else {
1825413bcc21SPierre Jolivet           mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1826413bcc21SPierre Jolivet         }
1827413bcc21SPierre Jolivet       }
1828413bcc21SPierre Jolivet       mumps->id.ICNTL(26) = -1;
1829413bcc21SPierre Jolivet     }
1830413bcc21SPierre Jolivet 
1831413bcc21SPierre Jolivet     /* copy MUMPS default control values from master to slaves. Although slaves do not call MUMPS, they may access these values in code.
1832413bcc21SPierre Jolivet        For example, ICNTL(9) is initialized to 1 by MUMPS and slaves check ICNTL(9) in MatSolve_MUMPS.
1833413bcc21SPierre Jolivet      */
1834413bcc21SPierre Jolivet     PetscCallMPI(MPI_Bcast(mumps->id.icntl, 40, MPI_INT, 0, mumps->omp_comm));
1835413bcc21SPierre Jolivet     PetscCallMPI(MPI_Bcast(mumps->id.cntl, 15, MPIU_REAL, 0, mumps->omp_comm));
1836413bcc21SPierre Jolivet 
1837413bcc21SPierre Jolivet     mumps->scat_rhs = NULL;
1838413bcc21SPierre Jolivet     mumps->scat_sol = NULL;
1839413bcc21SPierre Jolivet 
1840413bcc21SPierre Jolivet     /* set PETSc-MUMPS default options - override MUMPS default */
1841413bcc21SPierre Jolivet     mumps->id.ICNTL(3) = 0;
1842413bcc21SPierre Jolivet     mumps->id.ICNTL(4) = 0;
1843413bcc21SPierre Jolivet     if (mumps->petsc_size == 1) {
1844413bcc21SPierre Jolivet       mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */
1845413bcc21SPierre Jolivet       mumps->id.ICNTL(7)  = 7; /* automatic choice of ordering done by the package */
1846413bcc21SPierre Jolivet     } else {
1847413bcc21SPierre Jolivet       mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */
1848413bcc21SPierre Jolivet       mumps->id.ICNTL(21) = 1; /* distributed solution */
1849413bcc21SPierre Jolivet     }
1850413bcc21SPierre Jolivet   }
18519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_1", "ICNTL(1): output stream for error messages", "None", mumps->id.ICNTL(1), &icntl, &flg));
18529a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
18539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_2", "ICNTL(2): output stream for diagnostic printing, statistics, and warning", "None", mumps->id.ICNTL(2), &icntl, &flg));
18549a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
18559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_3", "ICNTL(3): output stream for global information, collected on the host", "None", mumps->id.ICNTL(3), &icntl, &flg));
18569a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
1857dcd589f8SShri Abhyankar 
18589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_4", "ICNTL(4): level of printing (0 to 4)", "None", mumps->id.ICNTL(4), &icntl, &flg));
18599a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
18609a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
18619a2535b5SHong Zhang 
18629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_6", "ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)", "None", mumps->id.ICNTL(6), &icntl, &flg));
18639a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
18649a2535b5SHong Zhang 
18659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_7", "ICNTL(7): computes a symmetric permutation in sequential analysis. 0=AMD, 2=AMF, 3=Scotch, 4=PORD, 5=Metis, 6=QAMD, and 7=auto(default)", "None", mumps->id.ICNTL(7), &icntl, &flg));
1866dcd589f8SShri Abhyankar   if (flg) {
1867aed4548fSBarry Smith     PetscCheck(icntl != 1 && icntl >= 0 && icntl <= 7, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Valid values are 0=AMD, 2=AMF, 3=Scotch, 4=PORD, 5=Metis, 6=QAMD, and 7=auto");
1868b53c1a7fSBarry Smith     mumps->id.ICNTL(7) = icntl;
1869dcd589f8SShri Abhyankar   }
1870e0b74bf9SHong Zhang 
18719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_8", "ICNTL(8): scaling strategy (-2 to 8 or 77)", "None", mumps->id.ICNTL(8), &mumps->id.ICNTL(8), NULL));
18729566063dSJacob Faibussowitsch   /* PetscCall(PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL)); handled by MatSolveTranspose_MUMPS() */
18739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_10", "ICNTL(10): max num of refinements", "None", mumps->id.ICNTL(10), &mumps->id.ICNTL(10), NULL));
18749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_11", "ICNTL(11): statistics related to an error analysis (via -ksp_view)", "None", mumps->id.ICNTL(11), &mumps->id.ICNTL(11), NULL));
18759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_12", "ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)", "None", mumps->id.ICNTL(12), &mumps->id.ICNTL(12), NULL));
18769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_13", "ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting", "None", mumps->id.ICNTL(13), &mumps->id.ICNTL(13), NULL));
18779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_14", "ICNTL(14): percentage increase in the estimated working space", "None", mumps->id.ICNTL(14), &mumps->id.ICNTL(14), NULL));
187845e3843bSPierre Jolivet   PetscCall(MatGetBlockSizes(A, &rbs, &cbs));
187945e3843bSPierre Jolivet   if (rbs == cbs && rbs > 1) mumps->id.ICNTL(15) = -rbs;
188045e3843bSPierre Jolivet   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_15", "ICNTL(15): compression of the input matrix resulting from a block format", "None", mumps->id.ICNTL(15), &mumps->id.ICNTL(15), &flg));
188145e3843bSPierre Jolivet   if (flg) {
188245e3843bSPierre Jolivet     PetscCheck(mumps->id.ICNTL(15) <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Positive -mat_mumps_icntl_15 not handled");
188345e3843bSPierre Jolivet     PetscCheck((-mumps->id.ICNTL(15) % cbs == 0) && (-mumps->id.ICNTL(15) % rbs == 0), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The opposite of -mat_mumps_icntl_15 must be a multiple of the column and row blocksizes");
188445e3843bSPierre Jolivet   }
18859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_19", "ICNTL(19): computes the Schur complement", "None", mumps->id.ICNTL(19), &mumps->id.ICNTL(19), NULL));
188659ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
18879566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F->schur));
18889566063dSJacob Faibussowitsch     PetscCall(MatMumpsResetSchur_Private(mumps));
188959ac8732SStefano Zampini   }
189025aac85cSJunchao Zhang 
189143f3b051SJunchao Zhang   /* Two MPICH Fortran MPI_IN_PLACE binding bugs prevented the use of 'mpich + mumps'. One happened with "mpi4py + mpich + mumps",
189243f3b051SJunchao Zhang      and was reported by Firedrake. See https://bitbucket.org/mpi4py/mpi4py/issues/162/mpi4py-initialization-breaks-fortran
189325aac85cSJunchao Zhang      and a petsc-maint mailing list thread with subject 'MUMPS segfaults in parallel because of ...'
189443f3b051SJunchao Zhang      This bug was fixed by https://github.com/pmodels/mpich/pull/4149. But the fix brought a new bug,
189543f3b051SJunchao Zhang      see https://github.com/pmodels/mpich/issues/5589. This bug was fixed by https://github.com/pmodels/mpich/pull/5590.
189643f3b051SJunchao Zhang      In short, we could not use distributed RHS with MPICH until v4.0b1.
189725aac85cSJunchao Zhang    */
189843f3b051SJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_LT(5, 3, 0) || (defined(PETSC_HAVE_MPICH_NUMVERSION) && (PETSC_HAVE_MPICH_NUMVERSION < 40000101))
189925aac85cSJunchao Zhang   mumps->ICNTL20 = 0; /* Centralized dense RHS*/
190043f3b051SJunchao Zhang #else
190143f3b051SJunchao Zhang   mumps->ICNTL20     = 10; /* Distributed dense RHS*/
190225aac85cSJunchao Zhang #endif
19039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_20", "ICNTL(20): give mumps centralized (0) or distributed (10) dense right-hand sides", "None", mumps->ICNTL20, &mumps->ICNTL20, &flg));
1904aed4548fSBarry Smith   PetscCheck(!flg || mumps->ICNTL20 == 10 || mumps->ICNTL20 == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "ICNTL(20)=%d is not supported by the PETSc/MUMPS interface. Allowed values are 0, 10", (int)mumps->ICNTL20);
190525aac85cSJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_LT(5, 3, 0)
1906aed4548fSBarry Smith   PetscCheck(!flg || mumps->ICNTL20 != 10, PETSC_COMM_SELF, PETSC_ERR_SUP, "ICNTL(20)=10 is not supported before MUMPS-5.3.0");
190725aac85cSJunchao Zhang #endif
19089566063dSJacob Faibussowitsch   /* PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL)); we only use distributed solution vector */
19099a2535b5SHong Zhang 
19109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_22", "ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)", "None", mumps->id.ICNTL(22), &mumps->id.ICNTL(22), NULL));
19119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_23", "ICNTL(23): max size of the working memory (MB) that can allocate per processor", "None", mumps->id.ICNTL(23), &mumps->id.ICNTL(23), NULL));
19129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_24", "ICNTL(24): detection of null pivot rows (0 or 1)", "None", mumps->id.ICNTL(24), &mumps->id.ICNTL(24), NULL));
19139371c9d4SSatish Balay   if (mumps->id.ICNTL(24)) { mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ }
1914d7ebd59bSHong Zhang 
19159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_25", "ICNTL(25): computes a solution of a deficient matrix and a null space basis", "None", mumps->id.ICNTL(25), &mumps->id.ICNTL(25), NULL));
19169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_26", "ICNTL(26): drives the solution phase if a Schur complement matrix", "None", mumps->id.ICNTL(26), &mumps->id.ICNTL(26), NULL));
19179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_27", "ICNTL(27): controls the blocking size for multiple right-hand sides", "None", mumps->id.ICNTL(27), &mumps->id.ICNTL(27), NULL));
19189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_28", "ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering", "None", mumps->id.ICNTL(28), &mumps->id.ICNTL(28), NULL));
19199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_29", "ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis", "None", mumps->id.ICNTL(29), &mumps->id.ICNTL(29), NULL));
19209566063dSJacob Faibussowitsch   /* PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL)); */ /* call MatMumpsGetInverse() directly */
19219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_31", "ICNTL(31): indicates which factors may be discarded during factorization", "None", mumps->id.ICNTL(31), &mumps->id.ICNTL(31), NULL));
19229566063dSJacob Faibussowitsch   /* PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL));  -- not supported by PETSc API */
19239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_33", "ICNTL(33): compute determinant", "None", mumps->id.ICNTL(33), &mumps->id.ICNTL(33), NULL));
19249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_35", "ICNTL(35): activates Block Low Rank (BLR) based factorization", "None", mumps->id.ICNTL(35), &mumps->id.ICNTL(35), NULL));
19259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_36", "ICNTL(36): choice of BLR factorization variant", "None", mumps->id.ICNTL(36), &mumps->id.ICNTL(36), NULL));
19269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_38", "ICNTL(38): estimated compression rate of LU factors with BLR", "None", mumps->id.ICNTL(38), &mumps->id.ICNTL(38), NULL));
1927dcd589f8SShri Abhyankar 
19289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_1", "CNTL(1): relative pivoting threshold", "None", mumps->id.CNTL(1), &mumps->id.CNTL(1), NULL));
19299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_2", "CNTL(2): stopping criterion of refinement", "None", mumps->id.CNTL(2), &mumps->id.CNTL(2), NULL));
19309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_3", "CNTL(3): absolute pivoting threshold", "None", mumps->id.CNTL(3), &mumps->id.CNTL(3), NULL));
19319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_4", "CNTL(4): value for static pivoting", "None", mumps->id.CNTL(4), &mumps->id.CNTL(4), NULL));
19329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_5", "CNTL(5): fixation for null pivots", "None", mumps->id.CNTL(5), &mumps->id.CNTL(5), NULL));
19339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_7", "CNTL(7): dropping parameter used during BLR", "None", mumps->id.CNTL(7), &mumps->id.CNTL(7), NULL));
1934e5bb22a1SHong Zhang 
19359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, sizeof(mumps->id.ooc_tmpdir), NULL));
1936b34f08ffSHong Zhang 
19379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-mat_mumps_view_info", "request INFO local to each processor", "", info, &ninfo, NULL));
1938b34f08ffSHong Zhang   if (ninfo) {
193908401ef6SPierre Jolivet     PetscCheck(ninfo <= 80, PETSC_COMM_SELF, PETSC_ERR_USER, "number of INFO %" PetscInt_FMT " must <= 80", ninfo);
19409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ninfo, &mumps->info));
1941b34f08ffSHong Zhang     mumps->ninfo = ninfo;
1942b34f08ffSHong Zhang     for (i = 0; i < ninfo; i++) {
1943aed4548fSBarry Smith       PetscCheck(info[i] >= 0 && info[i] <= 80, PETSC_COMM_SELF, PETSC_ERR_USER, "index of INFO %" PetscInt_FMT " must between 1 and 80", ninfo);
1944f7d195e4SLawrence Mitchell       mumps->info[i] = info[i];
1945b34f08ffSHong Zhang     }
1946b34f08ffSHong Zhang   }
1947d0609cedSBarry Smith   PetscOptionsEnd();
19483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1949dcd589f8SShri Abhyankar }
1950dcd589f8SShri Abhyankar 
1951d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F, Mat A, const MatFactorInfo *info, Mat_MUMPS *mumps)
1952d71ae5a4SJacob Faibussowitsch {
19535cd7cf9dSHong Zhang   PetscFunctionBegin;
19545cd7cf9dSHong Zhang   if (mumps->id.INFOG(1) < 0) {
19557a46b595SBarry Smith     PetscCheck(!A->erroriffailure, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in analysis phase: INFOG(1)=%d", mumps->id.INFOG(1));
19565cd7cf9dSHong Zhang     if (mumps->id.INFOG(1) == -6) {
19579566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1958603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
19595cd7cf9dSHong Zhang     } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
19609566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "problem of workspace, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1961603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1962dbf6bb8dSprj-     } else if (mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0) {
19639566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "Empty matrix\n"));
19645cd7cf9dSHong Zhang     } else {
19659566063dSJacob Faibussowitsch       PetscCall(PetscInfo(F, "Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
1966603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OTHER;
19675cd7cf9dSHong Zhang     }
19685cd7cf9dSHong Zhang   }
19693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19705cd7cf9dSHong Zhang }
19715cd7cf9dSHong Zhang 
1972d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
1973d71ae5a4SJacob Faibussowitsch {
1974e69c285eSBarry Smith   Mat_MUMPS     *mumps = (Mat_MUMPS *)F->data;
197567877ebaSShri Abhyankar   Vec            b;
197667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1977397b6df1SKris Buschelman 
1978397b6df1SKris Buschelman   PetscFunctionBegin;
1979d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
1980d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
19813ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1982d47f36abSHong Zhang   }
1983dcd589f8SShri Abhyankar 
19849a2535b5SHong Zhang   /* Set MUMPS options from the options database */
198526cc229bSBarry Smith   PetscCall(MatSetFromOptions_MUMPS(F, A));
1986dcd589f8SShri Abhyankar 
19879566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
19889566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX, mumps));
1989dcd589f8SShri Abhyankar 
199067877ebaSShri Abhyankar   /* analysis phase */
1991a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1992a5e57a09SHong Zhang   mumps->id.n   = M;
1993a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
199467877ebaSShri Abhyankar   case 0: /* centralized assembled matrix input */
1995a5e57a09SHong Zhang     if (!mumps->myid) {
1996a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
1997a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
1998a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
1999a6053eceSJunchao Zhang       if (mumps->id.ICNTL(6) > 1) mumps->id.a = (MumpsScalar *)mumps->val;
20004ac6704cSBarry Smith       if (r) {
20014ac6704cSBarry Smith         mumps->id.ICNTL(7) = 1;
2002a5e57a09SHong Zhang         if (!mumps->myid) {
2003e0b74bf9SHong Zhang           const PetscInt *idx;
2004a6053eceSJunchao Zhang           PetscInt        i;
20052205254eSKarl Rupp 
20069566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(M, &mumps->id.perm_in));
20079566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(r, &idx));
20089566063dSJacob Faibussowitsch           for (i = 0; i < M; i++) PetscCall(PetscMUMPSIntCast(idx[i] + 1, &(mumps->id.perm_in[i]))); /* perm_in[]: start from 1, not 0! */
20099566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(r, &idx));
2010e0b74bf9SHong Zhang         }
2011e0b74bf9SHong Zhang       }
201267877ebaSShri Abhyankar     }
201367877ebaSShri Abhyankar     break;
201467877ebaSShri Abhyankar   case 3: /* distributed assembled matrix input (size>1) */
2015a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2016a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2017a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2018a6053eceSJunchao Zhang     if (mumps->id.ICNTL(6) > 1) mumps->id.a_loc = (MumpsScalar *)mumps->val;
201925aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
20209566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
20219566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b, &mumps->scat_rhs, &mumps->b_seq));
20229566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
202325aac85cSJunchao Zhang     }
202467877ebaSShri Abhyankar     break;
202567877ebaSShri Abhyankar   }
20263ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
20279566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F, A, info, mumps));
202867877ebaSShri Abhyankar 
2029719d5645SBarry Smith   F->ops->lufactornumeric   = MatFactorNumeric_MUMPS;
2030dcd589f8SShri Abhyankar   F->ops->solve             = MatSolve_MUMPS;
203151d5961aSHong Zhang   F->ops->solvetranspose    = MatSolveTranspose_MUMPS;
20324e34a73bSHong Zhang   F->ops->matsolve          = MatMatSolve_MUMPS;
2033eb3ef3b2SHong Zhang   F->ops->mattransposesolve = MatMatTransposeSolve_MUMPS;
2034b18964edSHong Zhang   F->ops->matsolvetranspose = MatMatSolveTranspose_MUMPS;
2035d47f36abSHong Zhang 
2036d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
20373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2038b24902e0SBarry Smith }
2039b24902e0SBarry Smith 
2040450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
2041d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
2042d71ae5a4SJacob Faibussowitsch {
2043e69c285eSBarry Smith   Mat_MUMPS     *mumps = (Mat_MUMPS *)F->data;
204467877ebaSShri Abhyankar   Vec            b;
204567877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
2046450b117fSShri Abhyankar 
2047450b117fSShri Abhyankar   PetscFunctionBegin;
2048d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2049d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
20503ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2051d47f36abSHong Zhang   }
2052dcd589f8SShri Abhyankar 
20539a2535b5SHong Zhang   /* Set MUMPS options from the options database */
205426cc229bSBarry Smith   PetscCall(MatSetFromOptions_MUMPS(F, A));
2055dcd589f8SShri Abhyankar 
20569566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
20579566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX, mumps));
205867877ebaSShri Abhyankar 
205967877ebaSShri Abhyankar   /* analysis phase */
2060a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
2061a5e57a09SHong Zhang   mumps->id.n   = M;
2062a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
206367877ebaSShri Abhyankar   case 0: /* centralized assembled matrix input */
2064a5e57a09SHong Zhang     if (!mumps->myid) {
2065a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2066a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2067a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2068ad540459SPierre Jolivet       if (mumps->id.ICNTL(6) > 1) mumps->id.a = (MumpsScalar *)mumps->val;
206967877ebaSShri Abhyankar     }
207067877ebaSShri Abhyankar     break;
207167877ebaSShri Abhyankar   case 3: /* distributed assembled matrix input (size>1) */
2072a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2073a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2074a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2075ad540459SPierre Jolivet     if (mumps->id.ICNTL(6) > 1) mumps->id.a_loc = (MumpsScalar *)mumps->val;
207625aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
20779566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
20789566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b, &mumps->scat_rhs, &mumps->b_seq));
20799566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
208025aac85cSJunchao Zhang     }
208167877ebaSShri Abhyankar     break;
208267877ebaSShri Abhyankar   }
20833ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
20849566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F, A, info, mumps));
208567877ebaSShri Abhyankar 
2086450b117fSShri Abhyankar   F->ops->lufactornumeric   = MatFactorNumeric_MUMPS;
2087dcd589f8SShri Abhyankar   F->ops->solve             = MatSolve_MUMPS;
208851d5961aSHong Zhang   F->ops->solvetranspose    = MatSolveTranspose_MUMPS;
2089b18964edSHong Zhang   F->ops->matsolvetranspose = MatMatSolveTranspose_MUMPS;
2090d47f36abSHong Zhang 
2091d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
20923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2093450b117fSShri Abhyankar }
2094b24902e0SBarry Smith 
2095141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
2096d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F, Mat A, IS r, const MatFactorInfo *info)
2097d71ae5a4SJacob Faibussowitsch {
2098e69c285eSBarry Smith   Mat_MUMPS     *mumps = (Mat_MUMPS *)F->data;
209967877ebaSShri Abhyankar   Vec            b;
210067877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
2101397b6df1SKris Buschelman 
2102397b6df1SKris Buschelman   PetscFunctionBegin;
2103d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2104d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
21053ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2106d47f36abSHong Zhang   }
2107dcd589f8SShri Abhyankar 
21089a2535b5SHong Zhang   /* Set MUMPS options from the options database */
210926cc229bSBarry Smith   PetscCall(MatSetFromOptions_MUMPS(F, A));
2110dcd589f8SShri Abhyankar 
21119566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
21129566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX, mumps));
2113dcd589f8SShri Abhyankar 
211467877ebaSShri Abhyankar   /* analysis phase */
2115a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
2116a5e57a09SHong Zhang   mumps->id.n   = M;
2117a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
211867877ebaSShri Abhyankar   case 0: /* centralized assembled matrix input */
2119a5e57a09SHong Zhang     if (!mumps->myid) {
2120a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2121a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2122a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2123ad540459SPierre Jolivet       if (mumps->id.ICNTL(6) > 1) mumps->id.a = (MumpsScalar *)mumps->val;
212467877ebaSShri Abhyankar     }
212567877ebaSShri Abhyankar     break;
212667877ebaSShri Abhyankar   case 3: /* distributed assembled matrix input (size>1) */
2127a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2128a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2129a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2130ad540459SPierre Jolivet     if (mumps->id.ICNTL(6) > 1) mumps->id.a_loc = (MumpsScalar *)mumps->val;
213125aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
21329566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
21339566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b, &mumps->scat_rhs, &mumps->b_seq));
21349566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
213525aac85cSJunchao Zhang     }
213667877ebaSShri Abhyankar     break;
213767877ebaSShri Abhyankar   }
21383ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
21399566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F, A, info, mumps));
21405cd7cf9dSHong Zhang 
21412792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
2142dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
214351d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
21444e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
214523a5080aSHong Zhang   F->ops->mattransposesolve     = MatMatTransposeSolve_MUMPS;
2146b18964edSHong Zhang   F->ops->matsolvetranspose     = MatMatSolveTranspose_MUMPS;
21474e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
21480298fd71SBarry Smith   F->ops->getinertia = NULL;
21494e34a73bSHong Zhang #else
21504e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
2151db4efbfdSBarry Smith #endif
2152d47f36abSHong Zhang 
2153d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
21543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2155b24902e0SBarry Smith }
2156b24902e0SBarry Smith 
2157d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_MUMPS(Mat A, PetscViewer viewer)
2158d71ae5a4SJacob Faibussowitsch {
215964e6c443SBarry Smith   PetscBool         iascii;
216064e6c443SBarry Smith   PetscViewerFormat format;
2161e69c285eSBarry Smith   Mat_MUMPS        *mumps = (Mat_MUMPS *)A->data;
2162f6c57405SHong Zhang 
2163f6c57405SHong Zhang   PetscFunctionBegin;
216464e6c443SBarry Smith   /* check if matrix is mumps type */
21653ba16761SJacob Faibussowitsch   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(PETSC_SUCCESS);
216664e6c443SBarry Smith 
21679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
216864e6c443SBarry Smith   if (iascii) {
21699566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
21701511cd71SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
21719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "MUMPS run parameters:\n"));
21721511cd71SPierre Jolivet       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
21739566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  SYM (matrix type):                   %d\n", mumps->id.sym));
21749566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  PAR (host participation):            %d\n", mumps->id.par));
21759566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(1) (output for error):         %d\n", mumps->id.ICNTL(1)));
21769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(2) (output of diagnostic msg): %d\n", mumps->id.ICNTL(2)));
21779566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(3) (output for global info):   %d\n", mumps->id.ICNTL(3)));
21789566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(4) (level of printing):        %d\n", mumps->id.ICNTL(4)));
21799566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(5) (input mat struct):         %d\n", mumps->id.ICNTL(5)));
21809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(6) (matrix prescaling):        %d\n", mumps->id.ICNTL(6)));
21819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(7) (sequential matrix ordering):%d\n", mumps->id.ICNTL(7)));
21829566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(8) (scaling strategy):         %d\n", mumps->id.ICNTL(8)));
21839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(10) (max num of refinements):  %d\n", mumps->id.ICNTL(10)));
21849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(11) (error analysis):          %d\n", mumps->id.ICNTL(11)));
2185a5e57a09SHong Zhang         if (mumps->id.ICNTL(11) > 0) {
21869566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(4) (inf norm of input mat):        %g\n", mumps->id.RINFOG(4)));
21879566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(5) (inf norm of solution):         %g\n", mumps->id.RINFOG(5)));
21889566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(6) (inf norm of residual):         %g\n", mumps->id.RINFOG(6)));
21899566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n", mumps->id.RINFOG(7), mumps->id.RINFOG(8)));
21909566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(9) (error estimate):               %g\n", mumps->id.RINFOG(9)));
21919566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n", mumps->id.RINFOG(10), mumps->id.RINFOG(11)));
2192f6c57405SHong Zhang         }
21939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(12) (efficiency control):                         %d\n", mumps->id.ICNTL(12)));
21949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(13) (sequential factorization of the root node):  %d\n", mumps->id.ICNTL(13)));
21959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(14) (percentage of estimated workspace increase): %d\n", mumps->id.ICNTL(14)));
219645e3843bSPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(15) (compression of the input matrix):            %d\n", mumps->id.ICNTL(15)));
2197f6c57405SHong Zhang         /* ICNTL(15-17) not used */
21989566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(18) (input mat struct):                           %d\n", mumps->id.ICNTL(18)));
21999566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(19) (Schur complement info):                      %d\n", mumps->id.ICNTL(19)));
22009566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(20) (RHS sparse pattern):                         %d\n", mumps->id.ICNTL(20)));
22019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(21) (solution struct):                            %d\n", mumps->id.ICNTL(21)));
22029566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(22) (in-core/out-of-core facility):               %d\n", mumps->id.ICNTL(22)));
22039566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(23) (max size of memory can be allocated locally):%d\n", mumps->id.ICNTL(23)));
2204c0165424SHong Zhang 
22059566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(24) (detection of null pivot rows):               %d\n", mumps->id.ICNTL(24)));
22069566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(25) (computation of a null space basis):          %d\n", mumps->id.ICNTL(25)));
22079566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(26) (Schur options for RHS or solution):          %d\n", mumps->id.ICNTL(26)));
22089566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(27) (blocking size for multiple RHS):             %d\n", mumps->id.ICNTL(27)));
22099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(28) (use parallel or sequential ordering):        %d\n", mumps->id.ICNTL(28)));
22109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(29) (parallel ordering):                          %d\n", mumps->id.ICNTL(29)));
221142179a6aSHong Zhang 
22129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(30) (user-specified set of entries in inv(A)):    %d\n", mumps->id.ICNTL(30)));
22139566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(31) (factors is discarded in the solve phase):    %d\n", mumps->id.ICNTL(31)));
22149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(33) (compute determinant):                        %d\n", mumps->id.ICNTL(33)));
22159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(35) (activate BLR based factorization):           %d\n", mumps->id.ICNTL(35)));
22169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(36) (choice of BLR factorization variant):        %d\n", mumps->id.ICNTL(36)));
22179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(38) (estimated compression rate of LU factors):   %d\n", mumps->id.ICNTL(38)));
2218f6c57405SHong Zhang 
22199566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(1) (relative pivoting threshold):      %g\n", mumps->id.CNTL(1)));
22209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(2) (stopping criterion of refinement): %g\n", mumps->id.CNTL(2)));
22219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(3) (absolute pivoting threshold):      %g\n", mumps->id.CNTL(3)));
22229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(4) (value of static pivoting):         %g\n", mumps->id.CNTL(4)));
22239566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(5) (fixation for null pivots):         %g\n", mumps->id.CNTL(5)));
22249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(7) (dropping parameter for BLR):       %g\n", mumps->id.CNTL(7)));
2225f6c57405SHong Zhang 
2226a5b23f4aSJose E. Roman         /* information local to each processor */
22279566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis):\n"));
22289566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushSynchronized(viewer));
22299566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %g\n", mumps->myid, mumps->id.RINFO(1)));
22309566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
22319566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization):\n"));
22329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %g\n", mumps->myid, mumps->id.RINFO(2)));
22339566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
22349566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization):\n"));
22359566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %g\n", mumps->myid, mumps->id.RINFO(3)));
22369566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2237f6c57405SHong Zhang 
22389566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization):\n"));
22399566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(15)));
22409566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2241f6c57405SHong Zhang 
22429566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization):\n"));
22439566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(16)));
22449566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2245f6c57405SHong Zhang 
22469566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization):\n"));
22479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(23)));
22489566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2249b34f08ffSHong Zhang 
2250a0e18203SThibaut Appel         if (mumps->ninfo && mumps->ninfo <= 80) {
2251b34f08ffSHong Zhang           PetscInt i;
2252b34f08ffSHong Zhang           for (i = 0; i < mumps->ninfo; i++) {
22539566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(%" PetscInt_FMT "):\n", mumps->info[i]));
22549566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(mumps->info[i])));
22559566063dSJacob Faibussowitsch             PetscCall(PetscViewerFlush(viewer));
2256b34f08ffSHong Zhang           }
2257b34f08ffSHong Zhang         }
22589566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
22591511cd71SPierre Jolivet       } else PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all processes\n", ((PetscObject)A)->prefix ? ((PetscObject)A)->prefix : ""));
2260f6c57405SHong Zhang 
22611511cd71SPierre Jolivet       if (mumps->myid == 0) { /* information from the host */
22629566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFOG(1) (global estimated flops for the elimination after analysis): %g\n", mumps->id.RINFOG(1)));
22639566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFOG(2) (global estimated flops for the assembly after factorization): %g\n", mumps->id.RINFOG(2)));
22649566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFOG(3) (global estimated flops for the elimination after factorization): %g\n", mumps->id.RINFOG(3)));
22659566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n", mumps->id.RINFOG(12), mumps->id.RINFOG(13), mumps->id.INFOG(34)));
2266f6c57405SHong Zhang 
22679566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d\n", mumps->id.INFOG(3)));
22689566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d\n", mumps->id.INFOG(4)));
22699566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(5) (estimated maximum front size in the complete tree): %d\n", mumps->id.INFOG(5)));
22709566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(6) (number of nodes in the complete tree): %d\n", mumps->id.INFOG(6)));
22719566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(7) (ordering option effectively used after analysis): %d\n", mumps->id.INFOG(7)));
22729566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d\n", mumps->id.INFOG(8)));
22739566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d\n", mumps->id.INFOG(9)));
22749566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(10) (total integer space store the matrix factors after factorization): %d\n", mumps->id.INFOG(10)));
22759566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(11) (order of largest frontal matrix after factorization): %d\n", mumps->id.INFOG(11)));
22769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(12) (number of off-diagonal pivots): %d\n", mumps->id.INFOG(12)));
22779566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(13) (number of delayed pivots after factorization): %d\n", mumps->id.INFOG(13)));
22789566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(14) (number of memory compress after factorization): %d\n", mumps->id.INFOG(14)));
22799566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(15) (number of steps of iterative refinement after solution): %d\n", mumps->id.INFOG(15)));
22809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d\n", mumps->id.INFOG(16)));
22819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d\n", mumps->id.INFOG(17)));
22829566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d\n", mumps->id.INFOG(18)));
22839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d\n", mumps->id.INFOG(19)));
22849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(20) (estimated number of entries in the factors): %d\n", mumps->id.INFOG(20)));
22859566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d\n", mumps->id.INFOG(21)));
22869566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d\n", mumps->id.INFOG(22)));
22879566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d\n", mumps->id.INFOG(23)));
22889566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d\n", mumps->id.INFOG(24)));
22899566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d\n", mumps->id.INFOG(25)));
22909566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(28) (after factorization: number of null pivots encountered): %d\n", mumps->id.INFOG(28)));
22919566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n", mumps->id.INFOG(29)));
22929566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n", mumps->id.INFOG(30), mumps->id.INFOG(31)));
22939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(32) (after analysis: type of analysis done): %d\n", mumps->id.INFOG(32)));
22949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(33) (value used for ICNTL(8)): %d\n", mumps->id.INFOG(33)));
22959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(34) (exponent of the determinant if determinant is requested): %d\n", mumps->id.INFOG(34)));
22969566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(35) (after factorization: number of entries taking into account BLR factor compression - sum over all processors): %d\n", mumps->id.INFOG(35)));
22979566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(36) (after analysis: estimated size of all MUMPS internal data for running BLR in-core - value on the most memory consuming processor): %d\n", mumps->id.INFOG(36)));
22989566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(37) (after analysis: estimated size of all MUMPS internal data for running BLR in-core - sum over all processors): %d\n", mumps->id.INFOG(37)));
22999566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(38) (after analysis: estimated size of all MUMPS internal data for running BLR out-of-core - value on the most memory consuming processor): %d\n", mumps->id.INFOG(38)));
23009566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(39) (after analysis: estimated size of all MUMPS internal data for running BLR out-of-core - sum over all processors): %d\n", mumps->id.INFOG(39)));
2301f6c57405SHong Zhang       }
2302f6c57405SHong Zhang     }
2303cb828f0fSHong Zhang   }
23043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2305f6c57405SHong Zhang }
2306f6c57405SHong Zhang 
2307d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_MUMPS(Mat A, MatInfoType flag, MatInfo *info)
2308d71ae5a4SJacob Faibussowitsch {
2309e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)A->data;
231035bd34faSBarry Smith 
231135bd34faSBarry Smith   PetscFunctionBegin;
231235bd34faSBarry Smith   info->block_size        = 1.0;
2313cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
2314cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
231535bd34faSBarry Smith   info->nz_unneeded       = 0.0;
231635bd34faSBarry Smith   info->assemblies        = 0.0;
231735bd34faSBarry Smith   info->mallocs           = 0.0;
231835bd34faSBarry Smith   info->memory            = 0.0;
231935bd34faSBarry Smith   info->fill_ratio_given  = 0;
232035bd34faSBarry Smith   info->fill_ratio_needed = 0;
232135bd34faSBarry Smith   info->factor_mallocs    = 0;
23223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232335bd34faSBarry Smith }
232435bd34faSBarry Smith 
2325d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
2326d71ae5a4SJacob Faibussowitsch {
2327e69c285eSBarry Smith   Mat_MUMPS         *mumps = (Mat_MUMPS *)F->data;
2328a3d589ffSStefano Zampini   const PetscScalar *arr;
23298e7ba810SStefano Zampini   const PetscInt    *idxs;
23308e7ba810SStefano Zampini   PetscInt           size, i;
23316444a565SStefano Zampini 
23326444a565SStefano Zampini   PetscFunctionBegin;
23339566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &size));
2334b3cb21ddSStefano Zampini   /* Schur complement matrix */
23359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F->schur));
23369566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size, size, NULL, &F->schur));
23379566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(F->schur, &arr));
2338a3d589ffSStefano Zampini   mumps->id.schur      = (MumpsScalar *)arr;
2339a3d589ffSStefano Zampini   mumps->id.size_schur = size;
2340a3d589ffSStefano Zampini   mumps->id.schur_lld  = size;
23419566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(F->schur, &arr));
234248a46eb9SPierre Jolivet   if (mumps->sym == 1) PetscCall(MatSetOption(F->schur, MAT_SPD, PETSC_TRUE));
2343b3cb21ddSStefano Zampini 
2344b3cb21ddSStefano Zampini   /* MUMPS expects Fortran style indices */
23459566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.listvar_schur));
23469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &mumps->id.listvar_schur));
23479566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxs));
23489566063dSJacob Faibussowitsch   for (i = 0; i < size; i++) PetscCall(PetscMUMPSIntCast(idxs[i] + 1, &(mumps->id.listvar_schur[i])));
23499566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxs));
235059ac8732SStefano Zampini   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
2351b5fa320bSStefano Zampini   mumps->id.ICNTL(26) = -1;
23523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23536444a565SStefano Zampini }
235459ac8732SStefano Zampini 
2355d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F, Mat *S)
2356d71ae5a4SJacob Faibussowitsch {
23576444a565SStefano Zampini   Mat          St;
2358e69c285eSBarry Smith   Mat_MUMPS   *mumps = (Mat_MUMPS *)F->data;
23596444a565SStefano Zampini   PetscScalar *array;
23606444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX)
23618ac429a0SStefano Zampini   PetscScalar im = PetscSqrtScalar((PetscScalar)-1.0);
23626444a565SStefano Zampini #endif
23636444a565SStefano Zampini 
23646444a565SStefano Zampini   PetscFunctionBegin;
236508401ef6SPierre Jolivet   PetscCheck(mumps->id.ICNTL(19), PetscObjectComm((PetscObject)F), PETSC_ERR_ORDER, "Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
23669566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &St));
23679566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(St, PETSC_DECIDE, PETSC_DECIDE, mumps->id.size_schur, mumps->id.size_schur));
23689566063dSJacob Faibussowitsch   PetscCall(MatSetType(St, MATDENSE));
23699566063dSJacob Faibussowitsch   PetscCall(MatSetUp(St));
23709566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(St, &array));
237159ac8732SStefano Zampini   if (!mumps->sym) {                /* MUMPS always return a full matrix */
23726444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
23736444a565SStefano Zampini       PetscInt i, j, N = mumps->id.size_schur;
23746444a565SStefano Zampini       for (i = 0; i < N; i++) {
23756444a565SStefano Zampini         for (j = 0; j < N; j++) {
23766444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
23776444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j];
23786444a565SStefano Zampini #else
23796444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j].r + im * mumps->id.schur[i * N + j].i;
23806444a565SStefano Zampini #endif
23816444a565SStefano Zampini           array[j * N + i] = val;
23826444a565SStefano Zampini         }
23836444a565SStefano Zampini       }
23846444a565SStefano Zampini     } else { /* stored by columns */
23859566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array, mumps->id.schur, mumps->id.size_schur * mumps->id.size_schur));
23866444a565SStefano Zampini     }
23876444a565SStefano Zampini   } else {                          /* either full or lower-triangular (not packed) */
23886444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
23896444a565SStefano Zampini       PetscInt i, j, N = mumps->id.size_schur;
23906444a565SStefano Zampini       for (i = 0; i < N; i++) {
23916444a565SStefano Zampini         for (j = i; j < N; j++) {
23926444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
23936444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j];
23946444a565SStefano Zampini #else
23956444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j].r + im * mumps->id.schur[i * N + j].i;
23966444a565SStefano Zampini #endif
23976444a565SStefano Zampini           array[i * N + j] = val;
23986444a565SStefano Zampini           array[j * N + i] = val;
23996444a565SStefano Zampini         }
24006444a565SStefano Zampini       }
24016444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
24029566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array, mumps->id.schur, mumps->id.size_schur * mumps->id.size_schur));
24036444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
24046444a565SStefano Zampini       PetscInt i, j, N = mumps->id.size_schur;
24056444a565SStefano Zampini       for (i = 0; i < N; i++) {
24066444a565SStefano Zampini         for (j = 0; j < i + 1; j++) {
24076444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
24086444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j];
24096444a565SStefano Zampini #else
24106444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j].r + im * mumps->id.schur[i * N + j].i;
24116444a565SStefano Zampini #endif
24126444a565SStefano Zampini           array[i * N + j] = val;
24136444a565SStefano Zampini           array[j * N + i] = val;
24146444a565SStefano Zampini         }
24156444a565SStefano Zampini       }
24166444a565SStefano Zampini     }
24176444a565SStefano Zampini   }
24189566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(St, &array));
24196444a565SStefano Zampini   *S = St;
24203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24216444a565SStefano Zampini }
24226444a565SStefano Zampini 
2423d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F, PetscInt icntl, PetscInt ival)
2424d71ae5a4SJacob Faibussowitsch {
2425e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
24265ccb76cbSHong Zhang 
24275ccb76cbSHong Zhang   PetscFunctionBegin;
2428413bcc21SPierre Jolivet   if (mumps->id.job == JOB_NULL) {                                       /* need to cache icntl and ival since PetscMUMPS_c() has never been called */
2429413bcc21SPierre Jolivet     PetscInt i, nICNTL_pre = mumps->ICNTL_pre ? mumps->ICNTL_pre[0] : 0; /* number of already cached ICNTL */
24309371c9d4SSatish Balay     for (i = 0; i < nICNTL_pre; ++i)
24319371c9d4SSatish Balay       if (mumps->ICNTL_pre[1 + 2 * i] == icntl) break; /* is this ICNTL already cached? */
2432413bcc21SPierre Jolivet     if (i == nICNTL_pre) {                             /* not already cached */
2433413bcc21SPierre Jolivet       if (i > 0) PetscCall(PetscRealloc(sizeof(PetscMUMPSInt) * (2 * nICNTL_pre + 3), &mumps->ICNTL_pre));
2434413bcc21SPierre Jolivet       else PetscCall(PetscCalloc(sizeof(PetscMUMPSInt) * 3, &mumps->ICNTL_pre));
2435413bcc21SPierre Jolivet       mumps->ICNTL_pre[0]++;
2436413bcc21SPierre Jolivet     }
2437413bcc21SPierre Jolivet     mumps->ICNTL_pre[1 + 2 * i] = icntl;
2438413bcc21SPierre Jolivet     PetscCall(PetscMUMPSIntCast(ival, mumps->ICNTL_pre + 2 + 2 * i));
2439413bcc21SPierre Jolivet   } else PetscCall(PetscMUMPSIntCast(ival, &mumps->id.ICNTL(icntl)));
24403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24415ccb76cbSHong Zhang }
24425ccb76cbSHong Zhang 
2443d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F, PetscInt icntl, PetscInt *ival)
2444d71ae5a4SJacob Faibussowitsch {
2445e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2446bc6112feSHong Zhang 
2447bc6112feSHong Zhang   PetscFunctionBegin;
244836df9881Sjeremy theler   if (mumps->id.job == JOB_NULL) {
244936df9881Sjeremy theler     PetscInt i, nICNTL_pre = mumps->ICNTL_pre ? mumps->ICNTL_pre[0] : 0;
245036df9881Sjeremy theler     *ival = 0;
245136df9881Sjeremy theler     for (i = 0; i < nICNTL_pre; ++i) {
245236df9881Sjeremy theler       if (mumps->ICNTL_pre[1 + 2 * i] == icntl) *ival = mumps->ICNTL_pre[2 + 2 * i];
245336df9881Sjeremy theler     }
245436df9881Sjeremy theler   } else *ival = mumps->id.ICNTL(icntl);
24553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2456bc6112feSHong Zhang }
2457bc6112feSHong Zhang 
24585ccb76cbSHong Zhang /*@
24595ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
24605ccb76cbSHong Zhang 
2461c3339decSBarry Smith    Logically Collective
24625ccb76cbSHong Zhang 
24635ccb76cbSHong Zhang    Input Parameters:
246411a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
24655ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
24665ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
24675ccb76cbSHong Zhang 
24683c7db156SBarry Smith   Options Database Key:
2469147403d9SBarry Smith .   -mat_mumps_icntl_<icntl> <ival> - change the option numbered icntl to ival
24705ccb76cbSHong Zhang 
24715ccb76cbSHong Zhang    Level: beginner
24725ccb76cbSHong Zhang 
247396a0c994SBarry Smith    References:
2474606c0280SSatish Balay .  * - MUMPS Users' Guide
24755ccb76cbSHong Zhang 
2476*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
24775ccb76cbSHong Zhang @*/
2478d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsSetIcntl(Mat F, PetscInt icntl, PetscInt ival)
2479d71ae5a4SJacob Faibussowitsch {
24805ccb76cbSHong Zhang   PetscFunctionBegin;
24812989dfd4SHong Zhang   PetscValidType(F, 1);
248228b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
24835ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
24845ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F, ival, 3);
2485413bcc21SPierre Jolivet   PetscCheck(icntl >= 1 && icntl <= 38, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported ICNTL value %" PetscInt_FMT, icntl);
2486cac4c232SBarry Smith   PetscTryMethod(F, "MatMumpsSetIcntl_C", (Mat, PetscInt, PetscInt), (F, icntl, ival));
24873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24885ccb76cbSHong Zhang }
24895ccb76cbSHong Zhang 
2490a21f80fcSHong Zhang /*@
2491a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
2492a21f80fcSHong Zhang 
2493c3339decSBarry Smith    Logically Collective
2494a21f80fcSHong Zhang 
2495a21f80fcSHong Zhang    Input Parameters:
249611a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2497a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
2498a21f80fcSHong Zhang 
2499a21f80fcSHong Zhang   Output Parameter:
2500a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
2501a21f80fcSHong Zhang 
2502a21f80fcSHong Zhang    Level: beginner
2503a21f80fcSHong Zhang 
250496a0c994SBarry Smith    References:
2505606c0280SSatish Balay .  * - MUMPS Users' Guide
2506a21f80fcSHong Zhang 
2507*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
2508a21f80fcSHong Zhang @*/
2509d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetIcntl(Mat F, PetscInt icntl, PetscInt *ival)
2510d71ae5a4SJacob Faibussowitsch {
2511bc6112feSHong Zhang   PetscFunctionBegin;
25122989dfd4SHong Zhang   PetscValidType(F, 1);
251328b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2514bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
2515bc6112feSHong Zhang   PetscValidIntPointer(ival, 3);
2516413bcc21SPierre Jolivet   PetscCheck(icntl >= 1 && icntl <= 38, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported ICNTL value %" PetscInt_FMT, icntl);
2517cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetIcntl_C", (Mat, PetscInt, PetscInt *), (F, icntl, ival));
25183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2519bc6112feSHong Zhang }
2520bc6112feSHong Zhang 
2521d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F, PetscInt icntl, PetscReal val)
2522d71ae5a4SJacob Faibussowitsch {
2523e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
25248928b65cSHong Zhang 
25258928b65cSHong Zhang   PetscFunctionBegin;
2526413bcc21SPierre Jolivet   if (mumps->id.job == JOB_NULL) {
2527413bcc21SPierre Jolivet     PetscInt i, nCNTL_pre = mumps->CNTL_pre ? mumps->CNTL_pre[0] : 0;
25289371c9d4SSatish Balay     for (i = 0; i < nCNTL_pre; ++i)
25299371c9d4SSatish Balay       if (mumps->CNTL_pre[1 + 2 * i] == icntl) break;
2530413bcc21SPierre Jolivet     if (i == nCNTL_pre) {
2531413bcc21SPierre Jolivet       if (i > 0) PetscCall(PetscRealloc(sizeof(PetscReal) * (2 * nCNTL_pre + 3), &mumps->CNTL_pre));
2532413bcc21SPierre Jolivet       else PetscCall(PetscCalloc(sizeof(PetscReal) * 3, &mumps->CNTL_pre));
2533413bcc21SPierre Jolivet       mumps->CNTL_pre[0]++;
2534413bcc21SPierre Jolivet     }
2535413bcc21SPierre Jolivet     mumps->CNTL_pre[1 + 2 * i] = icntl;
2536413bcc21SPierre Jolivet     mumps->CNTL_pre[2 + 2 * i] = val;
2537413bcc21SPierre Jolivet   } else mumps->id.CNTL(icntl) = val;
25383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25398928b65cSHong Zhang }
25408928b65cSHong Zhang 
2541d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F, PetscInt icntl, PetscReal *val)
2542d71ae5a4SJacob Faibussowitsch {
2543e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2544bc6112feSHong Zhang 
2545bc6112feSHong Zhang   PetscFunctionBegin;
254636df9881Sjeremy theler   if (mumps->id.job == JOB_NULL) {
254736df9881Sjeremy theler     PetscInt i, nCNTL_pre = mumps->CNTL_pre ? mumps->CNTL_pre[0] : 0;
254836df9881Sjeremy theler     *val = 0.0;
254936df9881Sjeremy theler     for (i = 0; i < nCNTL_pre; ++i) {
255036df9881Sjeremy theler       if (mumps->CNTL_pre[1 + 2 * i] == icntl) *val = mumps->CNTL_pre[2 + 2 * i];
255136df9881Sjeremy theler     }
255236df9881Sjeremy theler   } else *val = mumps->id.CNTL(icntl);
25533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2554bc6112feSHong Zhang }
2555bc6112feSHong Zhang 
25568928b65cSHong Zhang /*@
25578928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
25588928b65cSHong Zhang 
2559c3339decSBarry Smith    Logically Collective
25608928b65cSHong Zhang 
25618928b65cSHong Zhang    Input Parameters:
256211a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
25638928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
25648928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
25658928b65cSHong Zhang 
25663c7db156SBarry Smith   Options Database Key:
2567147403d9SBarry Smith .   -mat_mumps_cntl_<icntl> <val>  - change the option numbered icntl to ival
25688928b65cSHong Zhang 
25698928b65cSHong Zhang    Level: beginner
25708928b65cSHong Zhang 
257196a0c994SBarry Smith    References:
2572606c0280SSatish Balay .  * - MUMPS Users' Guide
25738928b65cSHong Zhang 
2574*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
25758928b65cSHong Zhang @*/
2576d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsSetCntl(Mat F, PetscInt icntl, PetscReal val)
2577d71ae5a4SJacob Faibussowitsch {
25788928b65cSHong Zhang   PetscFunctionBegin;
25792989dfd4SHong Zhang   PetscValidType(F, 1);
258028b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
25818928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
2582bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F, val, 3);
2583413bcc21SPierre Jolivet   PetscCheck(icntl >= 1 && icntl <= 7, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported CNTL value %" PetscInt_FMT, icntl);
2584cac4c232SBarry Smith   PetscTryMethod(F, "MatMumpsSetCntl_C", (Mat, PetscInt, PetscReal), (F, icntl, val));
25853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25868928b65cSHong Zhang }
25878928b65cSHong Zhang 
2588a21f80fcSHong Zhang /*@
2589a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
2590a21f80fcSHong Zhang 
2591c3339decSBarry Smith    Logically Collective
2592a21f80fcSHong Zhang 
2593a21f80fcSHong Zhang    Input Parameters:
259411a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2595a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
2596a21f80fcSHong Zhang 
2597a21f80fcSHong Zhang   Output Parameter:
2598a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
2599a21f80fcSHong Zhang 
2600a21f80fcSHong Zhang    Level: beginner
2601a21f80fcSHong Zhang 
260296a0c994SBarry Smith    References:
2603606c0280SSatish Balay .  * - MUMPS Users' Guide
2604a21f80fcSHong Zhang 
2605*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
2606a21f80fcSHong Zhang @*/
2607d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetCntl(Mat F, PetscInt icntl, PetscReal *val)
2608d71ae5a4SJacob Faibussowitsch {
2609bc6112feSHong Zhang   PetscFunctionBegin;
26102989dfd4SHong Zhang   PetscValidType(F, 1);
261128b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2612bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
2613bc6112feSHong Zhang   PetscValidRealPointer(val, 3);
2614413bcc21SPierre Jolivet   PetscCheck(icntl >= 1 && icntl <= 7, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported CNTL value %" PetscInt_FMT, icntl);
2615cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetCntl_C", (Mat, PetscInt, PetscReal *), (F, icntl, val));
26163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2617bc6112feSHong Zhang }
2618bc6112feSHong Zhang 
2619d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F, PetscInt icntl, PetscInt *info)
2620d71ae5a4SJacob Faibussowitsch {
2621e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2622bc6112feSHong Zhang 
2623bc6112feSHong Zhang   PetscFunctionBegin;
2624bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
26253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2626bc6112feSHong Zhang }
2627bc6112feSHong Zhang 
2628d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F, PetscInt icntl, PetscInt *infog)
2629d71ae5a4SJacob Faibussowitsch {
2630e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2631bc6112feSHong Zhang 
2632bc6112feSHong Zhang   PetscFunctionBegin;
2633bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
26343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2635bc6112feSHong Zhang }
2636bc6112feSHong Zhang 
2637d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F, PetscInt icntl, PetscReal *rinfo)
2638d71ae5a4SJacob Faibussowitsch {
2639e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2640bc6112feSHong Zhang 
2641bc6112feSHong Zhang   PetscFunctionBegin;
2642bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
26433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2644bc6112feSHong Zhang }
2645bc6112feSHong Zhang 
2646d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F, PetscInt icntl, PetscReal *rinfog)
2647d71ae5a4SJacob Faibussowitsch {
2648e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2649bc6112feSHong Zhang 
2650bc6112feSHong Zhang   PetscFunctionBegin;
2651bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
26523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2653bc6112feSHong Zhang }
2654bc6112feSHong Zhang 
26555c0bae8cSAshish Patel PetscErrorCode MatMumpsGetNullPivots_MUMPS(Mat F, PetscInt *size, PetscInt **array)
26565c0bae8cSAshish Patel {
26575c0bae8cSAshish Patel   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
26585c0bae8cSAshish Patel 
26595c0bae8cSAshish Patel   PetscFunctionBegin;
26605c0bae8cSAshish Patel   PetscCheck(mumps->id.ICNTL(24) == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
26615c0bae8cSAshish Patel   *size  = 0;
26625c0bae8cSAshish Patel   *array = NULL;
26635c0bae8cSAshish Patel   if (!mumps->myid) {
26645c0bae8cSAshish Patel     *size = mumps->id.INFOG(28);
26655c0bae8cSAshish Patel     PetscCall(PetscMalloc1(*size, array));
26665c0bae8cSAshish Patel     for (int i = 0; i < *size; i++) (*array)[i] = mumps->id.pivnul_list[i] - 1;
26675c0bae8cSAshish Patel   }
26685c0bae8cSAshish Patel   PetscFunctionReturn(PETSC_SUCCESS);
26695c0bae8cSAshish Patel }
26705c0bae8cSAshish Patel 
2671d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F, Mat spRHS)
2672d71ae5a4SJacob Faibussowitsch {
26730e6b8875SHong Zhang   Mat          Bt = NULL, Btseq = NULL;
26740e6b8875SHong Zhang   PetscBool    flg;
2675bb599dfdSHong Zhang   Mat_MUMPS   *mumps = (Mat_MUMPS *)F->data;
2676bb599dfdSHong Zhang   PetscScalar *aa;
2677f410b75aSHong Zhang   PetscInt     spnr, *ia, *ja, M, nrhs;
2678bb599dfdSHong Zhang 
2679bb599dfdSHong Zhang   PetscFunctionBegin;
2680064a246eSJacob Faibussowitsch   PetscValidPointer(spRHS, 2);
2681013e2dc7SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)spRHS, MATTRANSPOSEVIRTUAL, &flg));
26820e6b8875SHong Zhang   if (flg) {
26839566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(spRHS, &Bt));
2684013e2dc7SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)spRHS), PETSC_ERR_ARG_WRONG, "Matrix spRHS must be type MATTRANSPOSEVIRTUAL matrix");
2685bb599dfdSHong Zhang 
26869566063dSJacob Faibussowitsch   PetscCall(MatMumpsSetIcntl(F, 30, 1));
2687bb599dfdSHong Zhang 
26882d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
26890e6b8875SHong Zhang     Mat_MPIAIJ *b = (Mat_MPIAIJ *)Bt->data;
26900e6b8875SHong Zhang     Btseq         = b->A;
26910e6b8875SHong Zhang   } else {
26920e6b8875SHong Zhang     Btseq = Bt;
26930e6b8875SHong Zhang   }
26940e6b8875SHong Zhang 
26959566063dSJacob Faibussowitsch   PetscCall(MatGetSize(spRHS, &M, &nrhs));
2696f410b75aSHong Zhang   mumps->id.nrhs = nrhs;
2697f410b75aSHong Zhang   mumps->id.lrhs = M;
2698f410b75aSHong Zhang   mumps->id.rhs  = NULL;
2699f410b75aSHong Zhang 
2700e3f2db6aSHong Zhang   if (!mumps->myid) {
27019566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArray(Btseq, &aa));
27029566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(Btseq, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
270328b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
27049566063dSJacob Faibussowitsch     PetscCall(PetscMUMPSIntCSRCast(mumps, spnr, ia, ja, &mumps->id.irhs_ptr, &mumps->id.irhs_sparse, &mumps->id.nz_rhs));
2705bb599dfdSHong Zhang     mumps->id.rhs_sparse = (MumpsScalar *)aa;
2706e3f2db6aSHong Zhang   } else {
2707e3f2db6aSHong Zhang     mumps->id.irhs_ptr    = NULL;
2708e3f2db6aSHong Zhang     mumps->id.irhs_sparse = NULL;
2709e3f2db6aSHong Zhang     mumps->id.nz_rhs      = 0;
2710e3f2db6aSHong Zhang     mumps->id.rhs_sparse  = NULL;
2711e3f2db6aSHong Zhang   }
2712bb599dfdSHong Zhang   mumps->id.ICNTL(20) = 1; /* rhs is sparse */
2713e3f2db6aSHong Zhang   mumps->id.ICNTL(21) = 0; /* solution is in assembled centralized format */
2714bb599dfdSHong Zhang 
2715bb599dfdSHong Zhang   /* solve phase */
2716bb599dfdSHong Zhang   mumps->id.job = JOB_SOLVE;
27173ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
2718049d1499SBarry Smith   PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MUMPS in solve phase: INFOG(1)=%d INFO(2)=%d", mumps->id.INFOG(1), mumps->id.INFO(2));
271914267174SHong Zhang 
2720e3f2db6aSHong Zhang   if (!mumps->myid) {
27219566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJRestoreArray(Btseq, &aa));
27229566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(Btseq, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
272328b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
2724e3f2db6aSHong Zhang   }
27253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2726bb599dfdSHong Zhang }
2727bb599dfdSHong Zhang 
2728bb599dfdSHong Zhang /*@
27292ef1f0ffSBarry Smith   MatMumpsGetInverse - Get user-specified set of entries in inverse of `A`
2730bb599dfdSHong Zhang 
2731c3339decSBarry Smith    Logically Collective
2732bb599dfdSHong Zhang 
273320f4b53cSBarry Smith    Input Parameter:
273420f4b53cSBarry Smith .  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2735bb599dfdSHong Zhang 
2736bb599dfdSHong Zhang   Output Parameter:
273720f4b53cSBarry Smith . spRHS - sequential sparse matrix in `MATTRANSPOSEVIRTUAL` format with requested entries of inverse of `A`
2738bb599dfdSHong Zhang 
2739bb599dfdSHong Zhang    Level: beginner
2740bb599dfdSHong Zhang 
2741bb599dfdSHong Zhang    References:
2742606c0280SSatish Balay .  * - MUMPS Users' Guide
2743bb599dfdSHong Zhang 
2744*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatCreateTranspose()`
2745bb599dfdSHong Zhang @*/
2746d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInverse(Mat F, Mat spRHS)
2747d71ae5a4SJacob Faibussowitsch {
2748bb599dfdSHong Zhang   PetscFunctionBegin;
2749bb599dfdSHong Zhang   PetscValidType(F, 1);
275028b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2751cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInverse_C", (Mat, Mat), (F, spRHS));
27523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2753bb599dfdSHong Zhang }
2754bb599dfdSHong Zhang 
2755d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInverseTranspose_MUMPS(Mat F, Mat spRHST)
2756d71ae5a4SJacob Faibussowitsch {
27570e6b8875SHong Zhang   Mat spRHS;
27580e6b8875SHong Zhang 
27590e6b8875SHong Zhang   PetscFunctionBegin;
27609566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(spRHST, &spRHS));
27619566063dSJacob Faibussowitsch   PetscCall(MatMumpsGetInverse_MUMPS(F, spRHS));
27629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&spRHS));
27633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27640e6b8875SHong Zhang }
27650e6b8875SHong Zhang 
27660e6b8875SHong Zhang /*@
27672ef1f0ffSBarry Smith   MatMumpsGetInverseTranspose - Get user-specified set of entries in inverse of matrix `A`^T
27680e6b8875SHong Zhang 
2769c3339decSBarry Smith    Logically Collective
27700e6b8875SHong Zhang 
277120f4b53cSBarry Smith    Input Parameter:
277220f4b53cSBarry Smith .  F - the factored matrix of A obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
27730e6b8875SHong Zhang 
27740e6b8875SHong Zhang   Output Parameter:
277520f4b53cSBarry Smith . spRHST - sequential sparse matrix in `MATAIJ` format containing the requested entries of inverse of `A`^T
27760e6b8875SHong Zhang 
27770e6b8875SHong Zhang    Level: beginner
27780e6b8875SHong Zhang 
27790e6b8875SHong Zhang    References:
2780606c0280SSatish Balay .  * - MUMPS Users' Guide
27810e6b8875SHong Zhang 
2782*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatCreateTranspose()`, `MatMumpsGetInverse()`
27830e6b8875SHong Zhang @*/
2784d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInverseTranspose(Mat F, Mat spRHST)
2785d71ae5a4SJacob Faibussowitsch {
27860e6b8875SHong Zhang   PetscBool flg;
27870e6b8875SHong Zhang 
27880e6b8875SHong Zhang   PetscFunctionBegin;
27890e6b8875SHong Zhang   PetscValidType(F, 1);
279028b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
27919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)spRHST, &flg, MATSEQAIJ, MATMPIAIJ, NULL));
279228b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)spRHST), PETSC_ERR_ARG_WRONG, "Matrix spRHST must be MATAIJ matrix");
27930e6b8875SHong Zhang 
2794cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInverseTranspose_C", (Mat, Mat), (F, spRHST));
27953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27960e6b8875SHong Zhang }
27970e6b8875SHong Zhang 
2798a21f80fcSHong Zhang /*@
2799a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
2800a21f80fcSHong Zhang 
2801c3339decSBarry Smith    Logically Collective
2802a21f80fcSHong Zhang 
2803a21f80fcSHong Zhang    Input Parameters:
280411a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2805a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
2806a21f80fcSHong Zhang 
2807a21f80fcSHong Zhang   Output Parameter:
2808a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
2809a21f80fcSHong Zhang 
2810a21f80fcSHong Zhang    Level: beginner
2811a21f80fcSHong Zhang 
281296a0c994SBarry Smith    References:
2813606c0280SSatish Balay .  * - MUMPS Users' Guide
2814a21f80fcSHong Zhang 
2815*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
2816a21f80fcSHong Zhang @*/
2817d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInfo(Mat F, PetscInt icntl, PetscInt *ival)
2818d71ae5a4SJacob Faibussowitsch {
2819bc6112feSHong Zhang   PetscFunctionBegin;
28202989dfd4SHong Zhang   PetscValidType(F, 1);
282128b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2822ca810319SHong Zhang   PetscValidIntPointer(ival, 3);
2823cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInfo_C", (Mat, PetscInt, PetscInt *), (F, icntl, ival));
28243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2825bc6112feSHong Zhang }
2826bc6112feSHong Zhang 
2827a21f80fcSHong Zhang /*@
2828a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2829a21f80fcSHong Zhang 
2830c3339decSBarry Smith    Logically Collective
2831a21f80fcSHong Zhang 
2832a21f80fcSHong Zhang    Input Parameters:
283311a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2834a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
2835a21f80fcSHong Zhang 
2836a21f80fcSHong Zhang   Output Parameter:
2837a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
2838a21f80fcSHong Zhang 
2839a21f80fcSHong Zhang    Level: beginner
2840a21f80fcSHong Zhang 
284196a0c994SBarry Smith    References:
2842606c0280SSatish Balay .  * - MUMPS Users' Guide
2843a21f80fcSHong Zhang 
2844*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
2845a21f80fcSHong Zhang @*/
2846d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInfog(Mat F, PetscInt icntl, PetscInt *ival)
2847d71ae5a4SJacob Faibussowitsch {
2848bc6112feSHong Zhang   PetscFunctionBegin;
28492989dfd4SHong Zhang   PetscValidType(F, 1);
285028b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2851ca810319SHong Zhang   PetscValidIntPointer(ival, 3);
2852cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInfog_C", (Mat, PetscInt, PetscInt *), (F, icntl, ival));
28533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2854bc6112feSHong Zhang }
2855bc6112feSHong Zhang 
2856a21f80fcSHong Zhang /*@
2857a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2858a21f80fcSHong Zhang 
2859c3339decSBarry Smith    Logically Collective
2860a21f80fcSHong Zhang 
2861a21f80fcSHong Zhang    Input Parameters:
286211a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2863a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
2864a21f80fcSHong Zhang 
2865a21f80fcSHong Zhang   Output Parameter:
2866a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
2867a21f80fcSHong Zhang 
2868a21f80fcSHong Zhang    Level: beginner
2869a21f80fcSHong Zhang 
287096a0c994SBarry Smith    References:
2871606c0280SSatish Balay .  * - MUMPS Users' Guide
2872a21f80fcSHong Zhang 
2873*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfog()`
2874a21f80fcSHong Zhang @*/
2875d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetRinfo(Mat F, PetscInt icntl, PetscReal *val)
2876d71ae5a4SJacob Faibussowitsch {
2877bc6112feSHong Zhang   PetscFunctionBegin;
28782989dfd4SHong Zhang   PetscValidType(F, 1);
287928b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2880bc6112feSHong Zhang   PetscValidRealPointer(val, 3);
2881cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetRinfo_C", (Mat, PetscInt, PetscReal *), (F, icntl, val));
28823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2883bc6112feSHong Zhang }
2884bc6112feSHong Zhang 
2885a21f80fcSHong Zhang /*@
2886a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2887a21f80fcSHong Zhang 
2888c3339decSBarry Smith    Logically Collective
2889a21f80fcSHong Zhang 
2890a21f80fcSHong Zhang    Input Parameters:
289111a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
2892a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
2893a21f80fcSHong Zhang 
2894a21f80fcSHong Zhang   Output Parameter:
2895a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
2896a21f80fcSHong Zhang 
2897a21f80fcSHong Zhang    Level: beginner
2898a21f80fcSHong Zhang 
289996a0c994SBarry Smith    References:
2900606c0280SSatish Balay .  * - MUMPS Users' Guide
2901a21f80fcSHong Zhang 
2902*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`
2903a21f80fcSHong Zhang @*/
2904d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetRinfog(Mat F, PetscInt icntl, PetscReal *val)
2905d71ae5a4SJacob Faibussowitsch {
2906bc6112feSHong Zhang   PetscFunctionBegin;
29072989dfd4SHong Zhang   PetscValidType(F, 1);
290828b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2909bc6112feSHong Zhang   PetscValidRealPointer(val, 3);
2910cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetRinfog_C", (Mat, PetscInt, PetscReal *), (F, icntl, val));
29113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2912bc6112feSHong Zhang }
2913bc6112feSHong Zhang 
29145c0bae8cSAshish Patel /*@
29155c0bae8cSAshish Patel   MatMumpsGetNullPivots - Get MUMPS parameter PIVNUL_LIST()
29165c0bae8cSAshish Patel 
29175c0bae8cSAshish Patel    Logically Collective
29185c0bae8cSAshish Patel 
29195c0bae8cSAshish Patel    Input Parameter:
29205c0bae8cSAshish Patel .  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-MUMPS interface
29215c0bae8cSAshish Patel 
29225c0bae8cSAshish Patel   Output Parameters:
29235c0bae8cSAshish Patel +  size - local size of the array. The size of the array is non-zero only on the host.
29245c0bae8cSAshish Patel -  array - array of rows with null pivot, these rows follow 0-based indexing. The array gets allocated within the function and the user is responsible
29255c0bae8cSAshish Patel            for freeing this array.
29265c0bae8cSAshish Patel 
29275c0bae8cSAshish Patel    Level: beginner
29285c0bae8cSAshish Patel 
29295c0bae8cSAshish Patel    References:
29305c0bae8cSAshish Patel .  * - MUMPS Users' Guide
29315c0bae8cSAshish Patel 
2932*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`
29335c0bae8cSAshish Patel @*/
29345c0bae8cSAshish Patel PetscErrorCode MatMumpsGetNullPivots(Mat F, PetscInt *size, PetscInt **array)
29355c0bae8cSAshish Patel {
29365c0bae8cSAshish Patel   PetscFunctionBegin;
29375c0bae8cSAshish Patel   PetscValidType(F, 1);
29385c0bae8cSAshish Patel   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
29395c0bae8cSAshish Patel   PetscValidIntPointer(size, 3);
29405c0bae8cSAshish Patel   PetscValidPointer(array, 4);
29415c0bae8cSAshish Patel   PetscUseMethod(F, "MatMumpsGetNullPivots_C", (Mat, PetscInt *, PetscInt **), (F, size, array));
29425c0bae8cSAshish Patel   PetscFunctionReturn(PETSC_SUCCESS);
29435c0bae8cSAshish Patel }
29445c0bae8cSAshish Patel 
294524b6179bSKris Buschelman /*MC
29462692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
294724b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
294824b6179bSKris Buschelman 
294911a5261eSBarry Smith   Works with `MATAIJ` and `MATSBAIJ` matrices
295024b6179bSKris Buschelman 
2951c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS
2952c2b89b5dSBarry Smith 
29532ef1f0ffSBarry Smith   Use ./configure --with-openmp --download-hwloc (or --with-hwloc) to enable running MUMPS in MPI+OpenMP hybrid mode and non-MUMPS in flat-MPI mode.
29542ef1f0ffSBarry Smith   See details below.
2955217d3b1eSJunchao Zhang 
29562ef1f0ffSBarry Smith   Use `-pc_type cholesky` or `lu` `-pc_factor_mat_solver_type mumps` to use this direct solver
2957c2b89b5dSBarry Smith 
295824b6179bSKris Buschelman   Options Database Keys:
29594422a9fcSPatrick Sanan +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
29604422a9fcSPatrick Sanan .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
29614422a9fcSPatrick Sanan .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
29624422a9fcSPatrick Sanan .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
29634422a9fcSPatrick Sanan .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2964b53c1a7fSBarry Smith .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis, 0=AMD, 2=AMF, 3=Scotch, 4=PORD, 5=Metis, 6=QAMD, and 7=auto
2965b53c1a7fSBarry Smith                         Use -pc_factor_mat_ordering_type <type> to have PETSc perform the ordering (sequential only)
29664422a9fcSPatrick Sanan .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
29674422a9fcSPatrick Sanan .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
29684422a9fcSPatrick Sanan .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
29694422a9fcSPatrick Sanan .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
29704422a9fcSPatrick Sanan .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
29714422a9fcSPatrick Sanan .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
297245e3843bSPierre Jolivet .  -mat_mumps_icntl_15  - ICNTL(15): compression of the input matrix resulting from a block format
29734422a9fcSPatrick Sanan .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
297425aac85cSJunchao Zhang .  -mat_mumps_icntl_20  - ICNTL(20): give MUMPS centralized (0) or distributed (10) dense RHS
29754422a9fcSPatrick Sanan .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
29764422a9fcSPatrick Sanan .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
29774422a9fcSPatrick Sanan .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
29784422a9fcSPatrick Sanan .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
29794422a9fcSPatrick Sanan .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
29804422a9fcSPatrick Sanan .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
29814422a9fcSPatrick Sanan .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
29824422a9fcSPatrick Sanan .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
29834422a9fcSPatrick Sanan .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
29844422a9fcSPatrick Sanan .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2985a0e18203SThibaut Appel .  -mat_mumps_icntl_35 - ICNTL(35): level of activation of BLR (Block Low-Rank) feature
2986a0e18203SThibaut Appel .  -mat_mumps_icntl_36 - ICNTL(36): controls the choice of BLR factorization variant
2987a0e18203SThibaut Appel .  -mat_mumps_icntl_38 - ICNTL(38): sets the estimated compression rate of LU factors with BLR
29884422a9fcSPatrick Sanan .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
29894422a9fcSPatrick Sanan .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
29904422a9fcSPatrick Sanan .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
29914422a9fcSPatrick Sanan .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2992217d3b1eSJunchao Zhang .  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2993a0e18203SThibaut Appel .  -mat_mumps_cntl_7 - CNTL(7): precision of the dropping parameter used during BLR factorization
2994217d3b1eSJunchao Zhang -  -mat_mumps_use_omp_threads [m] - run MUMPS in MPI+OpenMP hybrid mode as if omp_set_num_threads(m) is called before calling MUMPS.
2995217d3b1eSJunchao Zhang                                    Default might be the number of cores per CPU package (socket) as reported by hwloc and suggested by the MUMPS manual.
299624b6179bSKris Buschelman 
299724b6179bSKris Buschelman   Level: beginner
299824b6179bSKris Buschelman 
299995452b02SPatrick Sanan     Notes:
30002ef1f0ffSBarry Smith     MUMPS Cholesky does not handle (complex) Hermitian matrices (see User's Guide at https://mumps-solver.org/index.php?page=doc) so using it will
30012ef1f0ffSBarry Smith     error if the matrix is Hermitian.
300238548759SBarry Smith 
300326cc229bSBarry Smith     When used within a `KSP`/`PC` solve the options are prefixed with that of the `PC`. Otherwise one can set the options prefix by calling
300426cc229bSBarry Smith     `MatSetOptionsPrefixFactor()` on the matrix from which the factor was obtained or `MatSetOptionsPrefix()` on the factor matrix.
300526cc229bSBarry Smith 
30062ef1f0ffSBarry Smith     When a MUMPS factorization fails inside a KSP solve, for example with a `KSP_DIVERGED_PC_FAILED`, one can find the MUMPS information about
30072ef1f0ffSBarry Smith     the failure with
30082ef1f0ffSBarry Smith .vb
30092ef1f0ffSBarry Smith           KSPGetPC(ksp,&pc);
30102ef1f0ffSBarry Smith           PCFactorGetMatrix(pc,&mat);
30112ef1f0ffSBarry Smith           MatMumpsGetInfo(mat,....);
30122ef1f0ffSBarry Smith           MatMumpsGetInfog(mat,....); etc.
30132ef1f0ffSBarry Smith .ve
30142ef1f0ffSBarry Smith     Or run with `-ksp_error_if_not_converged` and the program will be stopped and the information printed in the error message.
30159fc87aa7SBarry Smith 
3016a5399872SJunchao Zhang     MUMPS provides 64-bit integer support in two build modes:
3017a5399872SJunchao Zhang       full 64-bit: here MUMPS is built with C preprocessing flag -DINTSIZE64 and Fortran compiler option -i8, -fdefault-integer-8 or equivalent, and
3018a5399872SJunchao Zhang       requires all dependent libraries MPI, ScaLAPACK, LAPACK and BLAS built the same way with 64-bit integers (for example ILP64 Intel MKL and MPI).
30198fcaa860SBarry Smith 
3020a5399872SJunchao Zhang       selective 64-bit: with the default MUMPS build, 64-bit integers have been introduced where needed. In compressed sparse row (CSR) storage of matrices,
3021a5399872SJunchao Zhang       MUMPS stores column indices in 32-bit, but row offsets in 64-bit, so you can have a huge number of non-zeros, but must have less than 2^31 rows and
3022a5399872SJunchao Zhang       columns. This can lead to significant memory and performance gains with respect to a full 64-bit integer MUMPS version. This requires a regular (32-bit
3023a5399872SJunchao Zhang       integer) build of all dependent libraries MPI, ScaLAPACK, LAPACK and BLAS.
3024a5399872SJunchao Zhang 
3025a5399872SJunchao Zhang     With --download-mumps=1, PETSc always build MUMPS in selective 64-bit mode, which can be used by both --with-64-bit-indices=0/1 variants of PETSc.
3026a5399872SJunchao Zhang 
3027a5399872SJunchao Zhang   Two modes to run MUMPS/PETSc with OpenMP
30282ef1f0ffSBarry Smith .vb
30292ef1f0ffSBarry Smith      Set OMP_NUM_THREADS and run with fewer MPI ranks than cores. For example, if you want to have 16 OpenMP
30302ef1f0ffSBarry Smith      threads per rank, then you may use "export OMP_NUM_THREADS=16 && mpirun -n 4 ./test".
30312ef1f0ffSBarry Smith .ve
30328fcaa860SBarry Smith 
30332ef1f0ffSBarry Smith .vb
30342ef1f0ffSBarry Smith      -mat_mumps_use_omp_threads [m] and run your code with as many MPI ranks as the number of cores. For example,
30352ef1f0ffSBarry Smith     if a compute node has 32 cores and you run on two nodes, you may use "mpirun -n 64 ./test -mat_mumps_use_omp_threads 16"
30362ef1f0ffSBarry Smith .ve
30378fcaa860SBarry Smith 
30388fcaa860SBarry Smith    To run MUMPS in MPI+OpenMP hybrid mode (i.e., enable multithreading in MUMPS), but still run the non-MUMPS part
30392ef1f0ffSBarry Smith    (i.e., PETSc part) of your code in the so-called flat-MPI (aka pure-MPI) mode, you need to configure PETSc with `--with-openmp` `--download-hwloc`
30402ef1f0ffSBarry Smith    (or `--with-hwloc`), and have an MPI that supports MPI-3.0's process shared memory (which is usually available). Since MUMPS calls BLAS
30418fcaa860SBarry Smith    libraries, to really get performance, you should have multithreaded BLAS libraries such as Intel MKL, AMD ACML, Cray libSci or OpenBLAS
30428fcaa860SBarry Smith    (PETSc will automatically try to utilized a threaded BLAS if --with-openmp is provided).
3043217d3b1eSJunchao Zhang 
30448fcaa860SBarry Smith    If you run your code through a job submission system, there are caveats in MPI rank mapping. We use MPI_Comm_split_type() to obtain MPI
3045217d3b1eSJunchao Zhang    processes on each compute node. Listing the processes in rank ascending order, we split processes on a node into consecutive groups of
3046217d3b1eSJunchao Zhang    size m and create a communicator called omp_comm for each group. Rank 0 in an omp_comm is called the master rank, and others in the omp_comm
3047217d3b1eSJunchao Zhang    are called slave ranks (or slaves). Only master ranks are seen to MUMPS and slaves are not. We will free CPUs assigned to slaves (might be set
3048217d3b1eSJunchao Zhang    by CPU binding policies in job scripts) and make the CPUs available to the master so that OMP threads spawned by MUMPS can run on the CPUs.
3049217d3b1eSJunchao Zhang    In a multi-socket compute node, MPI rank mapping is an issue. Still use the above example and suppose your compute node has two sockets,
3050217d3b1eSJunchao Zhang    if you interleave MPI ranks on the two sockets, in other words, even ranks are placed on socket 0, and odd ranks are on socket 1, and bind
3051217d3b1eSJunchao Zhang    MPI ranks to cores, then with -mat_mumps_use_omp_threads 16, a master rank (and threads it spawns) will use half cores in socket 0, and half
3052217d3b1eSJunchao Zhang    cores in socket 1, that definitely hurts locality. On the other hand, if you map MPI ranks consecutively on the two sockets, then the
3053217d3b1eSJunchao Zhang    problem will not happen. Therefore, when you use -mat_mumps_use_omp_threads, you need to keep an eye on your MPI rank mapping and CPU binding.
30548fcaa860SBarry Smith    For example, with the Slurm job scheduler, one can use srun --cpu-bind=verbose -m block:block to map consecutive MPI ranks to sockets and
3055217d3b1eSJunchao Zhang    examine the mapping result.
3056217d3b1eSJunchao Zhang 
305711a5261eSBarry Smith    PETSc does not control thread binding in MUMPS. So to get best performance, one still has to set `OMP_PROC_BIND` and `OMP_PLACES` in job scripts,
305811a5261eSBarry Smith    for example, export `OMP_PLACES`=threads and export `OMP_PROC_BIND`=spread. One does not need to export `OMP_NUM_THREADS`=m in job scripts as PETSc
305911a5261eSBarry Smith    calls `omp_set_num_threads`(m) internally before calling MUMPS.
3060217d3b1eSJunchao Zhang 
3061217d3b1eSJunchao Zhang    References:
3062606c0280SSatish Balay +  * - Heroux, Michael A., R. Brightwell, and Michael M. Wolf. "Bi-modal MPI and MPI+ threads computing on scalable multicore systems." IJHPCA (Submitted) (2011).
3063606c0280SSatish Balay -  * - Gutierrez, Samuel K., et al. "Accommodating Thread-Level Heterogeneity in Coupled Parallel Applications." Parallel and Distributed Processing Symposium (IPDPS), 2017 IEEE International. IEEE, 2017.
3064217d3b1eSJunchao Zhang 
3065*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`, `KSPGetPC()`, `PCFactorGetMatrix()`
306624b6179bSKris Buschelman M*/
306724b6179bSKris Buschelman 
3068d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_mumps(Mat A, MatSolverType *type)
3069d71ae5a4SJacob Faibussowitsch {
307035bd34faSBarry Smith   PetscFunctionBegin;
30712692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
30723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307335bd34faSBarry Smith }
307435bd34faSBarry Smith 
3075bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
3076d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_aij_mumps(Mat A, MatFactorType ftype, Mat *F)
3077d71ae5a4SJacob Faibussowitsch {
30782877fffaSHong Zhang   Mat         B;
30792877fffaSHong Zhang   Mat_MUMPS  *mumps;
3080ace3abfcSBarry Smith   PetscBool   isSeqAIJ;
30812c7c0729SBarry Smith   PetscMPIInt size;
30822877fffaSHong Zhang 
30832877fffaSHong Zhang   PetscFunctionBegin;
3084eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
3085b94d7dedSBarry Smith   PetscCheck(A->hermitian != PETSC_BOOL3_TRUE || A->symmetric == PETSC_BOOL3_TRUE || ftype != MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_SUP, "Hermitian CHOLESKY Factor is not supported");
3086eb1ec7c1SStefano Zampini #endif
30872877fffaSHong Zhang   /* Create the factorization matrix */
30889566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeqAIJ));
30899566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
30909566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
30919566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps", &((PetscObject)B)->type_name));
30929566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
30932877fffaSHong Zhang 
30944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
30952205254eSKarl Rupp 
30962877fffaSHong Zhang   B->ops->view    = MatView_MUMPS;
309735bd34faSBarry Smith   B->ops->getinfo = MatGetInfo_MUMPS;
30982205254eSKarl Rupp 
30999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
31009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
31019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
31029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
31039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
31049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
31059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
31069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
31079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
31089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
31099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
31105c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
31119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
31129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
31136444a565SStefano Zampini 
3114450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
3115450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
3116d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
3117bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
3118bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
31199566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
3120746480a1SHong Zhang     mumps->sym = 0;
3121dcd589f8SShri Abhyankar   } else {
312267877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
3123450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
3124bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
3125bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
31269566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
312759ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
312859ac8732SStefano Zampini     mumps->sym = 2;
312959ac8732SStefano Zampini #else
3130b94d7dedSBarry Smith     if (A->spd == PETSC_BOOL3_TRUE) mumps->sym = 1;
31316fdc2a6dSBarry Smith     else mumps->sym = 2;
313259ac8732SStefano Zampini #endif
3133450b117fSShri Abhyankar   }
31342877fffaSHong Zhang 
313500c67f3bSHong Zhang   /* set solvertype */
31369566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
31379566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
31389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
31392c7c0729SBarry Smith   if (size == 1) {
31404ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3141f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
31422c7c0729SBarry Smith   }
31432877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
3144e69c285eSBarry Smith   B->data         = (void *)mumps;
31452205254eSKarl Rupp 
31462877fffaSHong Zhang   *F               = B;
3147413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3148413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3149413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3150d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
31513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31522877fffaSHong Zhang }
31532877fffaSHong Zhang 
3154bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
3155d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A, MatFactorType ftype, Mat *F)
3156d71ae5a4SJacob Faibussowitsch {
31572877fffaSHong Zhang   Mat         B;
31582877fffaSHong Zhang   Mat_MUMPS  *mumps;
3159ace3abfcSBarry Smith   PetscBool   isSeqSBAIJ;
31602c7c0729SBarry Smith   PetscMPIInt size;
31612877fffaSHong Zhang 
31622877fffaSHong Zhang   PetscFunctionBegin;
3163eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
3164b94d7dedSBarry Smith   PetscCheck(A->hermitian != PETSC_BOOL3_TRUE || A->symmetric == PETSC_BOOL3_TRUE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Hermitian CHOLESKY Factor is not supported");
3165eb1ec7c1SStefano Zampini #endif
31669566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
31679566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
31689566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps", &((PetscObject)B)->type_name));
31699566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3170e69c285eSBarry Smith 
31714dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
31729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSeqSBAIJ));
3173bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
317416ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
3175dcd589f8SShri Abhyankar   } else {
3176bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
3177bccb9932SShri Abhyankar   }
3178bccb9932SShri Abhyankar 
317967877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
3180bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
3181722b6324SPierre Jolivet   B->ops->getinfo                = MatGetInfo_MUMPS;
31822205254eSKarl Rupp 
31839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
31849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
31859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
31869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
31879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
31889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
31899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
31909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
31919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
31929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
31939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
31945c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
31959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
31969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
31972205254eSKarl Rupp 
3198f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
319959ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
320059ac8732SStefano Zampini   mumps->sym = 2;
320159ac8732SStefano Zampini #else
3202b94d7dedSBarry Smith   if (A->spd == PETSC_BOOL3_TRUE) mumps->sym = 1;
32036fdc2a6dSBarry Smith   else mumps->sym = 2;
320459ac8732SStefano Zampini #endif
3205a214ac2aSShri Abhyankar 
320600c67f3bSHong Zhang   /* set solvertype */
32079566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
32089566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
32099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
32102c7c0729SBarry Smith   if (size == 1) {
32114ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3212f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
32132c7c0729SBarry Smith   }
32149566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
3215f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
3216e69c285eSBarry Smith   B->data         = (void *)mumps;
32172205254eSKarl Rupp 
32182877fffaSHong Zhang   *F               = B;
3219413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3220413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3221413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3222d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
32233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32242877fffaSHong Zhang }
322597969023SHong Zhang 
3226d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_baij_mumps(Mat A, MatFactorType ftype, Mat *F)
3227d71ae5a4SJacob Faibussowitsch {
322867877ebaSShri Abhyankar   Mat         B;
322967877ebaSShri Abhyankar   Mat_MUMPS  *mumps;
3230ace3abfcSBarry Smith   PetscBool   isSeqBAIJ;
32312c7c0729SBarry Smith   PetscMPIInt size;
323267877ebaSShri Abhyankar 
323367877ebaSShri Abhyankar   PetscFunctionBegin;
323467877ebaSShri Abhyankar   /* Create the factorization matrix */
32359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &isSeqBAIJ));
32369566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
32379566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
32389566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps", &((PetscObject)B)->type_name));
32399566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3240450b117fSShri Abhyankar 
32414dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
3242450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
3243450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
3244450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
3245bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
3246bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
3247746480a1SHong Zhang     mumps->sym = 0;
32489566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
3249546078acSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead");
3250bccb9932SShri Abhyankar 
3251450b117fSShri Abhyankar   B->ops->view    = MatView_MUMPS;
3252722b6324SPierre Jolivet   B->ops->getinfo = MatGetInfo_MUMPS;
32532205254eSKarl Rupp 
32549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
32559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
32569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
32579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
32589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
32599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
32609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
32619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
32629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
32639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
32649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
32655c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
32669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
32679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
3268450b117fSShri Abhyankar 
326900c67f3bSHong Zhang   /* set solvertype */
32709566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
32719566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
32729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
32732c7c0729SBarry Smith   if (size == 1) {
32744ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3275f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
32762c7c0729SBarry Smith   }
32777ee00b23SStefano Zampini   B->ops->destroy = MatDestroy_MUMPS;
32787ee00b23SStefano Zampini   B->data         = (void *)mumps;
32797ee00b23SStefano Zampini 
32807ee00b23SStefano Zampini   *F               = B;
3281413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3282413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3283413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3284d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
32853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32867ee00b23SStefano Zampini }
32877ee00b23SStefano Zampini 
32887ee00b23SStefano Zampini /* MatGetFactor for Seq and MPI SELL matrices */
3289d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_sell_mumps(Mat A, MatFactorType ftype, Mat *F)
3290d71ae5a4SJacob Faibussowitsch {
32917ee00b23SStefano Zampini   Mat         B;
32927ee00b23SStefano Zampini   Mat_MUMPS  *mumps;
32937ee00b23SStefano Zampini   PetscBool   isSeqSELL;
32942c7c0729SBarry Smith   PetscMPIInt size;
32957ee00b23SStefano Zampini 
32967ee00b23SStefano Zampini   PetscFunctionBegin;
32977ee00b23SStefano Zampini   /* Create the factorization matrix */
32989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSELL, &isSeqSELL));
32999566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
33009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
33019566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps", &((PetscObject)B)->type_name));
33029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
33037ee00b23SStefano Zampini 
33044dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
33057ee00b23SStefano Zampini 
33067ee00b23SStefano Zampini   B->ops->view    = MatView_MUMPS;
33077ee00b23SStefano Zampini   B->ops->getinfo = MatGetInfo_MUMPS;
33087ee00b23SStefano Zampini 
33099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
33109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
33119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
33129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
33139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
33149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
33159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
33169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
33179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
33189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
33199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
33205c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
33217ee00b23SStefano Zampini 
33227ee00b23SStefano Zampini   if (ftype == MAT_FACTOR_LU) {
33237ee00b23SStefano Zampini     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
33247ee00b23SStefano Zampini     B->factortype            = MAT_FACTOR_LU;
33257ee00b23SStefano Zampini     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
33267ee00b23SStefano Zampini     else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "To be implemented");
33277ee00b23SStefano Zampini     mumps->sym = 0;
33289566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
33297ee00b23SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "To be implemented");
33307ee00b23SStefano Zampini 
33317ee00b23SStefano Zampini   /* set solvertype */
33329566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
33339566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
33349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
33352c7c0729SBarry Smith   if (size == 1) {
33364ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization  */
3337f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
33382c7c0729SBarry Smith   }
3339450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
3340e69c285eSBarry Smith   B->data         = (void *)mumps;
33412205254eSKarl Rupp 
3342450b117fSShri Abhyankar   *F               = B;
3343413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3344413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3345413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3346d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
33473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3348450b117fSShri Abhyankar }
334942c9c57cSBarry Smith 
3350d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
3351d71ae5a4SJacob Faibussowitsch {
335242c9c57cSBarry Smith   PetscFunctionBegin;
33539566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
33549566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
33559566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIBAIJ, MAT_FACTOR_LU, MatGetFactor_baij_mumps));
33569566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_baij_mumps));
33579566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPISBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_sbaij_mumps));
33589566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
33599566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
33609566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQBAIJ, MAT_FACTOR_LU, MatGetFactor_baij_mumps));
33619566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_baij_mumps));
33629566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_sbaij_mumps));
33639566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQSELL, MAT_FACTOR_LU, MatGetFactor_sell_mumps));
33643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
336542c9c57cSBarry Smith }
3366