xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 50ea2040ea806f6bc490422fc45783ec8be71b35)
1397b6df1SKris Buschelman /*
2c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
3397b6df1SKris Buschelman */
467602552SJunchao Zhang #include <petscpkg_version.h>
59d0448ceSStefano Zampini #include <petscsf.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 
109261f6e4SBarry Smith #define MUMPS_MANUALS "(see users manual https://mumps-solver.org/index.php?page=doc \"Error and warning diagnostics\")"
119261f6e4SBarry Smith 
12397b6df1SKris Buschelman EXTERN_C_BEGIN
13397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
142907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
152907cef9SHong Zhang     #include <cmumps_c.h>
162907cef9SHong Zhang   #else
17c6db04a5SJed Brown     #include <zmumps_c.h>
182907cef9SHong Zhang   #endif
192907cef9SHong Zhang #else
202907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
212907cef9SHong Zhang     #include <smumps_c.h>
22397b6df1SKris Buschelman   #else
23c6db04a5SJed Brown     #include <dmumps_c.h>
24397b6df1SKris Buschelman   #endif
252907cef9SHong Zhang #endif
26397b6df1SKris Buschelman EXTERN_C_END
27397b6df1SKris Buschelman #define JOB_INIT         -1
28413bcc21SPierre Jolivet #define JOB_NULL         0
293d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
303d472b54SHong Zhang #define JOB_FACTNUMERIC  2
313d472b54SHong Zhang #define JOB_SOLVE        3
32397b6df1SKris Buschelman #define JOB_END          -2
333d472b54SHong Zhang 
342907cef9SHong Zhang /* calls to MUMPS */
352907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
362907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
373ab56b82SJunchao Zhang     #define MUMPS_c cmumps_c
382907cef9SHong Zhang   #else
393ab56b82SJunchao Zhang     #define MUMPS_c zmumps_c
402907cef9SHong Zhang   #endif
412907cef9SHong Zhang #else
422907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
433ab56b82SJunchao Zhang     #define MUMPS_c smumps_c
442907cef9SHong Zhang   #else
453ab56b82SJunchao Zhang     #define MUMPS_c dmumps_c
462907cef9SHong Zhang   #endif
472907cef9SHong Zhang #endif
482907cef9SHong Zhang 
49a6053eceSJunchao Zhang /* MUMPS uses MUMPS_INT for nonzero indices such as irn/jcn, irn_loc/jcn_loc and uses int64_t for
50a6053eceSJunchao Zhang    number of nonzeros such as nnz, nnz_loc. We typedef MUMPS_INT to PetscMUMPSInt to follow the
51a6053eceSJunchao Zhang    naming convention in PetscMPIInt, PetscBLASInt etc.
52a6053eceSJunchao Zhang */
53a6053eceSJunchao Zhang typedef MUMPS_INT PetscMUMPSInt;
54a6053eceSJunchao Zhang 
5567602552SJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_GE(5, 3, 0)
5667602552SJunchao 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 */
57a6053eceSJunchao Zhang     #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
5867602552SJunchao Zhang   #endif
59a6053eceSJunchao Zhang #else
6067602552SJunchao Zhang   #if defined(INTSIZE64) /* INTSIZE64 is a command line macro one used to build MUMPS in full 64-bit mode */
6167602552SJunchao Zhang     #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
6267602552SJunchao Zhang   #endif
6367602552SJunchao Zhang #endif
6467602552SJunchao Zhang 
65a6053eceSJunchao Zhang #define MPIU_MUMPSINT       MPI_INT
66a6053eceSJunchao Zhang #define PETSC_MUMPS_INT_MAX 2147483647
67a6053eceSJunchao Zhang #define PETSC_MUMPS_INT_MIN -2147483648
68a6053eceSJunchao Zhang 
69a6053eceSJunchao Zhang /* Cast PetscInt to PetscMUMPSInt. Usually there is no overflow since <a> is row/col indices or some small integers*/
706497c311SBarry Smith static inline PetscErrorCode PetscMUMPSIntCast(PetscCount a, PetscMUMPSInt *b)
71d71ae5a4SJacob Faibussowitsch {
72a6053eceSJunchao Zhang   PetscFunctionBegin;
73ece88022SPierre Jolivet #if PetscDefined(USE_64BIT_INDICES)
742c71b3e2SJacob Faibussowitsch   PetscAssert(a <= PETSC_MUMPS_INT_MAX && a >= PETSC_MUMPS_INT_MIN, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "PetscInt too long for PetscMUMPSInt");
75ece88022SPierre Jolivet #endif
7657508eceSPierre Jolivet   *b = (PetscMUMPSInt)a;
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78a6053eceSJunchao Zhang }
79a6053eceSJunchao Zhang 
80a6053eceSJunchao Zhang /* Put these utility routines here since they are only used in this file */
81ce78bad3SBarry Smith 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)
82d71ae5a4SJacob Faibussowitsch {
83a6053eceSJunchao Zhang   PetscInt  myval;
84a6053eceSJunchao Zhang   PetscBool myset;
854d86920dSPierre Jolivet 
86a6053eceSJunchao Zhang   PetscFunctionBegin;
87a6053eceSJunchao Zhang   /* PetscInt's size should be always >= PetscMUMPSInt's. It is safe to call PetscOptionsInt_Private to read a PetscMUMPSInt */
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, (PetscInt)currentvalue, &myval, &myset, lb, ub));
899566063dSJacob Faibussowitsch   if (myset) PetscCall(PetscMUMPSIntCast(myval, value));
90a6053eceSJunchao Zhang   if (set) *set = myset;
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92a6053eceSJunchao Zhang }
93a6053eceSJunchao 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)
94a6053eceSJunchao Zhang 
95217d3b1eSJunchao 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 */
963ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
973ab56b82SJunchao Zhang   #define PetscMUMPS_c(mumps) \
983ab56b82SJunchao Zhang     do { \
993ab56b82SJunchao Zhang       if (mumps->use_petsc_omp_support) { \
1003ab56b82SJunchao Zhang         if (mumps->is_omp_master) { \
1019566063dSJacob Faibussowitsch           PetscCall(PetscOmpCtrlOmpRegionOnMasterBegin(mumps->omp_ctrl)); \
10214ffdc6fSStefano Zampini           PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
10314ffdc6fSStefano Zampini           PetscStackCallExternalVoid(PetscStringize(MUMPS_c), MUMPS_c(&mumps->id)); \
10414ffdc6fSStefano Zampini           PetscCall(PetscFPTrapPop()); \
1059566063dSJacob Faibussowitsch           PetscCall(PetscOmpCtrlOmpRegionOnMasterEnd(mumps->omp_ctrl)); \
1063ab56b82SJunchao Zhang         } \
1079566063dSJacob Faibussowitsch         PetscCall(PetscOmpCtrlBarrier(mumps->omp_ctrl)); \
108c3714a1dSJunchao Zhang         /* Global info is same on all processes so we Bcast it within omp_comm. Local info is specific      \
109c3714a1dSJunchao Zhang          to processes, so we only Bcast info[1], an error code and leave others (since they do not have   \
110c3714a1dSJunchao Zhang          an easy translation between omp_comm and petsc_comm). See MUMPS-5.1.2 manual p82.                   \
111c3714a1dSJunchao Zhang          omp_comm is a small shared memory communicator, hence doing multiple Bcast as shown below is OK. \
112c3714a1dSJunchao Zhang       */ \
113338d3105SPierre Jolivet         PetscCallMPI(MPI_Bcast(mumps->id.infog, PETSC_STATIC_ARRAY_LENGTH(mumps->id.infog), MPIU_MUMPSINT, 0, mumps->omp_comm)); \
114338d3105SPierre Jolivet         PetscCallMPI(MPI_Bcast(mumps->id.rinfog, PETSC_STATIC_ARRAY_LENGTH(mumps->id.rinfog), MPIU_REAL, 0, mumps->omp_comm)); \
115338d3105SPierre Jolivet         PetscCallMPI(MPI_Bcast(mumps->id.info, PETSC_STATIC_ARRAY_LENGTH(mumps->id.info), MPIU_MUMPSINT, 0, mumps->omp_comm)); \
116338d3105SPierre Jolivet         PetscCallMPI(MPI_Bcast(mumps->id.rinfo, PETSC_STATIC_ARRAY_LENGTH(mumps->id.rinfo), MPIU_REAL, 0, mumps->omp_comm)); \
1173ab56b82SJunchao Zhang       } else { \
11814ffdc6fSStefano Zampini         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
11914ffdc6fSStefano Zampini         PetscStackCallExternalVoid(PetscStringize(MUMPS_c), MUMPS_c(&mumps->id)); \
12014ffdc6fSStefano Zampini         PetscCall(PetscFPTrapPop()); \
1213ab56b82SJunchao Zhang       } \
1223ab56b82SJunchao Zhang     } while (0)
1233ab56b82SJunchao Zhang #else
1243ab56b82SJunchao Zhang   #define PetscMUMPS_c(mumps) \
125d71ae5a4SJacob Faibussowitsch     do { \
12614ffdc6fSStefano Zampini       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
12714ffdc6fSStefano Zampini       PetscStackCallExternalVoid(PetscStringize(MUMPS_c), MUMPS_c(&mumps->id)); \
12814ffdc6fSStefano Zampini       PetscCall(PetscFPTrapPop()); \
129d71ae5a4SJacob Faibussowitsch     } while (0)
1303ab56b82SJunchao Zhang #endif
1313ab56b82SJunchao Zhang 
132940cd9d6SSatish Balay /* declare MumpsScalar */
133940cd9d6SSatish Balay #if defined(PETSC_USE_COMPLEX)
134940cd9d6SSatish Balay   #if defined(PETSC_USE_REAL_SINGLE)
135940cd9d6SSatish Balay     #define MumpsScalar mumps_complex
136940cd9d6SSatish Balay   #else
137940cd9d6SSatish Balay     #define MumpsScalar mumps_double_complex
138940cd9d6SSatish Balay   #endif
139940cd9d6SSatish Balay #else
140940cd9d6SSatish Balay   #define MumpsScalar PetscScalar
141940cd9d6SSatish Balay #endif
1423d472b54SHong Zhang 
143397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
144397b6df1SKris Buschelman #define ICNTL(I)  icntl[(I) - 1]
145397b6df1SKris Buschelman #define CNTL(I)   cntl[(I) - 1]
146397b6df1SKris Buschelman #define INFOG(I)  infog[(I) - 1]
147a7aca84bSHong Zhang #define INFO(I)   info[(I) - 1]
148397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I) - 1]
149adc1d99fSHong Zhang #define RINFO(I)  rinfo[(I) - 1]
150397b6df1SKris Buschelman 
151a6053eceSJunchao Zhang typedef struct Mat_MUMPS Mat_MUMPS;
152a6053eceSJunchao Zhang struct Mat_MUMPS {
153397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
1542907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
1552907cef9SHong Zhang   CMUMPS_STRUC_C id;
1562907cef9SHong Zhang   #else
157397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
1582907cef9SHong Zhang   #endif
1592907cef9SHong Zhang #else
1602907cef9SHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
1612907cef9SHong Zhang   SMUMPS_STRUC_C id;
162397b6df1SKris Buschelman   #else
163397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
164397b6df1SKris Buschelman   #endif
1652907cef9SHong Zhang #endif
1662907cef9SHong Zhang 
167397b6df1SKris Buschelman   MatStructure   matstruc;
1682d4298aeSJunchao Zhang   PetscMPIInt    myid, petsc_size;
169a6053eceSJunchao Zhang   PetscMUMPSInt *irn, *jcn;       /* the (i,j,v) triplets passed to mumps. */
170a6053eceSJunchao 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. */
1716497c311SBarry Smith   PetscCount     nnz;             /* number of nonzeros. The type is called selective 64-bit in mumps */
172a6053eceSJunchao Zhang   PetscMUMPSInt  sym;
1732d4298aeSJunchao Zhang   MPI_Comm       mumps_comm;
174413bcc21SPierre Jolivet   PetscMUMPSInt *ICNTL_pre;
175413bcc21SPierre Jolivet   PetscReal     *CNTL_pre;
176a6053eceSJunchao Zhang   PetscMUMPSInt  ICNTL9_pre;         /* check if ICNTL(9) is changed from previous MatSolve */
177801fbe65SHong Zhang   VecScatter     scat_rhs, scat_sol; /* used by MatSolve() */
17825aac85cSJunchao Zhang   PetscMUMPSInt  ICNTL20;            /* use centralized (0) or distributed (10) dense RHS */
17967602552SJunchao Zhang   PetscMUMPSInt  lrhs_loc, nloc_rhs, *irhs_loc;
18067602552SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
18167602552SJunchao Zhang   PetscInt    *rhs_nrow, max_nrhs;
18267602552SJunchao Zhang   PetscMPIInt *rhs_recvcounts, *rhs_disps;
18367602552SJunchao Zhang   PetscScalar *rhs_loc, *rhs_recvbuf;
18467602552SJunchao Zhang #endif
185801fbe65SHong Zhang   Vec            b_seq, x_seq;
186a6053eceSJunchao Zhang   PetscInt       ninfo, *info; /* which INFO to display */
187b5fa320bSStefano Zampini   PetscInt       sizeredrhs;
18859ac8732SStefano Zampini   PetscScalar   *schur_sol;
18959ac8732SStefano Zampini   PetscInt       schur_sizesol;
190a6053eceSJunchao Zhang   PetscMUMPSInt *ia_alloc, *ja_alloc; /* work arrays used for the CSR struct for sparse rhs */
1916497c311SBarry Smith   PetscCount     cur_ilen, cur_jlen;  /* current len of ia_alloc[], ja_alloc[] */
192a6053eceSJunchao Zhang   PetscErrorCode (*ConvertToTriples)(Mat, PetscInt, MatReuse, Mat_MUMPS *);
1932205254eSKarl Rupp 
1949d0448ceSStefano Zampini   /* Support for MATNEST */
1959d0448ceSStefano Zampini   PetscErrorCode (**nest_convert_to_triples)(Mat, PetscInt, MatReuse, Mat_MUMPS *);
1966497c311SBarry Smith   PetscCount  *nest_vals_start;
1979d0448ceSStefano Zampini   PetscScalar *nest_vals;
1989d0448ceSStefano Zampini 
199a6053eceSJunchao Zhang   /* stuff used by petsc/mumps OpenMP support*/
2003ab56b82SJunchao Zhang   PetscBool    use_petsc_omp_support;
201da81f932SPierre Jolivet   PetscOmpCtrl omp_ctrl;             /* an OpenMP controller that blocked processes will release their CPU (MPI_Barrier does not have this guarantee) */
2023ab56b82SJunchao Zhang   MPI_Comm     petsc_comm, omp_comm; /* petsc_comm is petsc matrix's comm */
2036497c311SBarry Smith   PetscCount  *recvcount;            /* a collection of nnz on omp_master */
204a6053eceSJunchao Zhang   PetscMPIInt  tag, omp_comm_size;
2053ab56b82SJunchao Zhang   PetscBool    is_omp_master; /* is this rank the master of omp_comm */
206a6053eceSJunchao Zhang   MPI_Request *reqs;
207a6053eceSJunchao Zhang };
2083ab56b82SJunchao Zhang 
209a6053eceSJunchao Zhang /* Cast a 1-based CSR represented by (nrow, ia, ja) of type PetscInt to a CSR of type PetscMUMPSInt.
210a6053eceSJunchao Zhang    Here, nrow is number of rows, ia[] is row pointer and ja[] is column indices.
211a6053eceSJunchao Zhang  */
212d2a308c1SPierre Jolivet static PetscErrorCode PetscMUMPSIntCSRCast(PETSC_UNUSED Mat_MUMPS *mumps, PetscInt nrow, PetscInt *ia, PetscInt *ja, PetscMUMPSInt **ia_mumps, PetscMUMPSInt **ja_mumps, PetscMUMPSInt *nnz_mumps)
213d71ae5a4SJacob Faibussowitsch {
2146497c311SBarry Smith   PetscInt nnz = ia[nrow] - 1; /* mumps uses 1-based indices. Uses PetscInt instead of PetscCount since mumps only uses PetscMUMPSInt for rhs */
215f0c56d0fSKris Buschelman 
216a6053eceSJunchao Zhang   PetscFunctionBegin;
217a6053eceSJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
218a6053eceSJunchao Zhang   {
219a6053eceSJunchao Zhang     PetscInt i;
220a6053eceSJunchao Zhang     if (nrow + 1 > mumps->cur_ilen) { /* realloc ia_alloc/ja_alloc to fit ia/ja */
2219566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->ia_alloc));
2229566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrow + 1, &mumps->ia_alloc));
223a6053eceSJunchao Zhang       mumps->cur_ilen = nrow + 1;
224a6053eceSJunchao Zhang     }
225a6053eceSJunchao Zhang     if (nnz > mumps->cur_jlen) {
2269566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->ja_alloc));
2279566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nnz, &mumps->ja_alloc));
228a6053eceSJunchao Zhang       mumps->cur_jlen = nnz;
229a6053eceSJunchao Zhang     }
230f4f49eeaSPierre Jolivet     for (i = 0; i < nrow + 1; i++) PetscCall(PetscMUMPSIntCast(ia[i], &mumps->ia_alloc[i]));
231f4f49eeaSPierre Jolivet     for (i = 0; i < nnz; i++) PetscCall(PetscMUMPSIntCast(ja[i], &mumps->ja_alloc[i]));
232a6053eceSJunchao Zhang     *ia_mumps = mumps->ia_alloc;
233a6053eceSJunchao Zhang     *ja_mumps = mumps->ja_alloc;
234a6053eceSJunchao Zhang   }
235a6053eceSJunchao Zhang #else
236a6053eceSJunchao Zhang   *ia_mumps = ia;
237a6053eceSJunchao Zhang   *ja_mumps = ja;
238a6053eceSJunchao Zhang #endif
2399566063dSJacob Faibussowitsch   PetscCall(PetscMUMPSIntCast(nnz, nnz_mumps));
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241a6053eceSJunchao Zhang }
242b24902e0SBarry Smith 
243d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS *mumps)
244d71ae5a4SJacob Faibussowitsch {
245b5fa320bSStefano Zampini   PetscFunctionBegin;
2469566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.listvar_schur));
2479566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.redrhs));
2489566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->schur_sol));
24959ac8732SStefano Zampini   mumps->id.size_schur = 0;
250b3cb21ddSStefano Zampini   mumps->id.schur_lld  = 0;
25159ac8732SStefano Zampini   mumps->id.ICNTL(19)  = 0;
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25359ac8732SStefano Zampini }
25459ac8732SStefano Zampini 
255b3cb21ddSStefano Zampini /* solve with rhs in mumps->id.redrhs and return in the same location */
256d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
257d71ae5a4SJacob Faibussowitsch {
258b3cb21ddSStefano Zampini   Mat_MUMPS           *mumps = (Mat_MUMPS *)F->data;
259b3cb21ddSStefano Zampini   Mat                  S, B, X;
260b3cb21ddSStefano Zampini   MatFactorSchurStatus schurstatus;
261b3cb21ddSStefano Zampini   PetscInt             sizesol;
26259ac8732SStefano Zampini 
26359ac8732SStefano Zampini   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(MatFactorFactorizeSchurComplement(F));
2659566063dSJacob Faibussowitsch   PetscCall(MatFactorGetSchurComplement(F, &S, &schurstatus));
2669566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mumps->id.size_schur, mumps->id.nrhs, (PetscScalar *)mumps->id.redrhs, &B));
2679566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, ((PetscObject)S)->type_name));
268a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2699566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(B, S->boundtocpu));
270a3d589ffSStefano Zampini #endif
271b3cb21ddSStefano Zampini   switch (schurstatus) {
272d71ae5a4SJacob Faibussowitsch   case MAT_FACTOR_SCHUR_FACTORED:
273d71ae5a4SJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mumps->id.size_schur, mumps->id.nrhs, (PetscScalar *)mumps->id.redrhs, &X));
274d71ae5a4SJacob Faibussowitsch     PetscCall(MatSetType(X, ((PetscObject)S)->type_name));
275a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2769566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(X, S->boundtocpu));
277a3d589ffSStefano Zampini #endif
278b3cb21ddSStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
2799566063dSJacob Faibussowitsch       PetscCall(MatMatSolveTranspose(S, B, X));
28059ac8732SStefano Zampini     } else {
2819566063dSJacob Faibussowitsch       PetscCall(MatMatSolve(S, B, X));
28259ac8732SStefano Zampini     }
283b3cb21ddSStefano Zampini     break;
284b3cb21ddSStefano Zampini   case MAT_FACTOR_SCHUR_INVERTED:
285b3cb21ddSStefano Zampini     sizesol = mumps->id.nrhs * mumps->id.size_schur;
28659ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
2879566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->schur_sol));
2889566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sizesol, &mumps->schur_sol));
28959ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
290b5fa320bSStefano Zampini     }
2919566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mumps->id.size_schur, mumps->id.nrhs, mumps->schur_sol, &X));
2929566063dSJacob Faibussowitsch     PetscCall(MatSetType(X, ((PetscObject)S)->type_name));
293a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2949566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(X, S->boundtocpu));
295a3d589ffSStefano Zampini #endif
2969566063dSJacob Faibussowitsch     PetscCall(MatProductCreateWithMat(S, B, NULL, X));
29759ac8732SStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
2989566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(X, MATPRODUCT_AtB));
299b5fa320bSStefano Zampini     } else {
3009566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(X, MATPRODUCT_AB));
301b5fa320bSStefano Zampini     }
3029566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(X));
3039566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(X));
3049566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(X));
3054417c5e8SHong Zhang 
3069566063dSJacob Faibussowitsch     PetscCall(MatCopy(X, B, SAME_NONZERO_PATTERN));
307b3cb21ddSStefano Zampini     break;
308d71ae5a4SJacob Faibussowitsch   default:
309d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status);
31059ac8732SStefano Zampini   }
3119566063dSJacob Faibussowitsch   PetscCall(MatFactorRestoreSchurComplement(F, &S, schurstatus));
3129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
3139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
315b5fa320bSStefano Zampini }
316b5fa320bSStefano Zampini 
317d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
318d71ae5a4SJacob Faibussowitsch {
319b3cb21ddSStefano Zampini   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
320b5fa320bSStefano Zampini 
321b5fa320bSStefano Zampini   PetscFunctionBegin;
322b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
3233ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
324b5fa320bSStefano Zampini   }
325b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
326b5fa320bSStefano Zampini     PetscInt sizeredrhs = mumps->id.nrhs * mumps->id.size_schur;
327b5fa320bSStefano Zampini     /* allocate MUMPS internal array to store reduced right-hand sides */
328b5fa320bSStefano Zampini     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
3299566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->id.redrhs));
330b5fa320bSStefano Zampini       mumps->id.lredrhs = mumps->id.size_schur;
3319566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mumps->id.nrhs * mumps->id.lredrhs, &mumps->id.redrhs));
332b5fa320bSStefano Zampini       mumps->sizeredrhs = mumps->id.nrhs * mumps->id.lredrhs;
333b5fa320bSStefano Zampini     }
334b5fa320bSStefano Zampini   } else { /* prepare for the expansion step */
335b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
3369566063dSJacob Faibussowitsch     PetscCall(MatMumpsSolveSchur_Private(F));
337b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
3383ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
3399261f6e4SBarry Smith     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in solve: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
340b5fa320bSStefano Zampini     /* restore defaults */
341b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
342d3d598ffSStefano Zampini     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
343d3d598ffSStefano Zampini     if (mumps->id.nrhs > 1) {
3449566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->id.redrhs));
345d3d598ffSStefano Zampini       mumps->id.lredrhs = 0;
346d3d598ffSStefano Zampini       mumps->sizeredrhs = 0;
347d3d598ffSStefano Zampini     }
348b5fa320bSStefano Zampini   }
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
350b5fa320bSStefano Zampini }
351b5fa320bSStefano Zampini 
352397b6df1SKris Buschelman /*
353d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
354d341cd04SHong Zhang 
355397b6df1SKris Buschelman   input:
35675480915SPierre Jolivet     A       - matrix in aij,baij or sbaij format
357397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
358bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
359bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
360397b6df1SKris Buschelman   output:
361397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
362397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
363eb9baa12SBarry Smith 
364eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
3657ee00b23SStefano Zampini   freed with PetscFree(mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
366eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
367eb9baa12SBarry Smith 
368397b6df1SKris Buschelman  */
36916ebf90aSShri Abhyankar 
37066976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
371d71ae5a4SJacob Faibussowitsch {
372a3d589ffSStefano Zampini   const PetscScalar *av;
373185f6596SHong Zhang   const PetscInt    *ai, *aj, *ajj, M = A->rmap->n;
3746497c311SBarry Smith   PetscCount         nz, rnz, k;
375a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
37616ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa = (Mat_SeqAIJ *)A->data;
377397b6df1SKris Buschelman 
378397b6df1SKris Buschelman   PetscFunctionBegin;
3799566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(A, &av));
380bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
3812205254eSKarl Rupp     nz = aa->nz;
3822205254eSKarl Rupp     ai = aa->i;
3832205254eSKarl Rupp     aj = aa->j;
3849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
3856497c311SBarry Smith     for (PetscCount i = k = 0; i < M; i++) {
38616ebf90aSShri Abhyankar       rnz = ai[i + 1] - ai[i];
38767877ebaSShri Abhyankar       ajj = aj + ai[i];
3886497c311SBarry Smith       for (PetscCount j = 0; j < rnz; j++) {
3899566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(i + shift, &row[k]));
3909566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[k]));
391a6053eceSJunchao Zhang         k++;
39216ebf90aSShri Abhyankar       }
39316ebf90aSShri Abhyankar     }
39450c845baSStefano Zampini     mumps->val = (PetscScalar *)av;
395a6053eceSJunchao Zhang     mumps->irn = row;
396a6053eceSJunchao Zhang     mumps->jcn = col;
397a6053eceSJunchao Zhang     mumps->nnz = nz;
398127cd276SPierre Jolivet   } else if (mumps->nest_vals) PetscCall(PetscArraycpy(mumps->val, av, aa->nz)); /* MatConvertToTriples_nest_xaij() allocates mumps->val outside of MatConvertToTriples_seqaij_seqaij(), so one needs to copy the memory */
399127cd276SPierre Jolivet   else mumps->val = (PetscScalar *)av;                                           /* in the default case, mumps->val is never allocated, one just needs to update the mumps->val pointer */
4009566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40216ebf90aSShri Abhyankar }
403397b6df1SKris Buschelman 
40466976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
405d71ae5a4SJacob Faibussowitsch {
4066497c311SBarry Smith   PetscCount     nz, i, j, k, r;
4077ee00b23SStefano Zampini   Mat_SeqSELL   *a = (Mat_SeqSELL *)A->data;
408a6053eceSJunchao Zhang   PetscMUMPSInt *row, *col;
4097ee00b23SStefano Zampini 
4107ee00b23SStefano Zampini   PetscFunctionBegin;
4117ee00b23SStefano Zampini   nz = a->sliidx[a->totalslices];
41250c845baSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
4139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
414a6053eceSJunchao Zhang     for (i = k = 0; i < a->totalslices; i++) {
41548a46eb9SPierre 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++]));
4167ee00b23SStefano Zampini     }
4179566063dSJacob Faibussowitsch     for (i = 0; i < nz; i++) PetscCall(PetscMUMPSIntCast(a->colidx[i] + shift, &col[i]));
418a6053eceSJunchao Zhang     mumps->irn = row;
419a6053eceSJunchao Zhang     mumps->jcn = col;
420a6053eceSJunchao Zhang     mumps->nnz = nz;
42150c845baSStefano Zampini     mumps->val = a->val;
422127cd276SPierre Jolivet   } else if (mumps->nest_vals) PetscCall(PetscArraycpy(mumps->val, a->val, nz)); /* MatConvertToTriples_nest_xaij() allocates mumps->val outside of MatConvertToTriples_seqsell_seqaij(), so one needs to copy the memory */
423127cd276SPierre Jolivet   else mumps->val = a->val;                                                      /* in the default case, mumps->val is never allocated, one just needs to update the mumps->val pointer */
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4257ee00b23SStefano Zampini }
4267ee00b23SStefano Zampini 
42766976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
428d71ae5a4SJacob Faibussowitsch {
42967877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa = (Mat_SeqBAIJ *)A->data;
43033d57670SJed Brown   const PetscInt *ai, *aj, *ajj, bs2 = aa->bs2;
4316497c311SBarry Smith   PetscCount      M, nz = bs2 * aa->nz, idx = 0, rnz, i, j, k, m;
432a6053eceSJunchao Zhang   PetscInt        bs;
433a6053eceSJunchao Zhang   PetscMUMPSInt  *row, *col;
43467877ebaSShri Abhyankar 
43567877ebaSShri Abhyankar   PetscFunctionBegin;
43650c845baSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
4379566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A, &bs));
43833d57670SJed Brown     M  = A->rmap->N / bs;
4399371c9d4SSatish Balay     ai = aa->i;
4409371c9d4SSatish Balay     aj = aa->j;
4419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
44267877ebaSShri Abhyankar     for (i = 0; i < M; i++) {
44367877ebaSShri Abhyankar       ajj = aj + ai[i];
44467877ebaSShri Abhyankar       rnz = ai[i + 1] - ai[i];
44567877ebaSShri Abhyankar       for (k = 0; k < rnz; k++) {
44667877ebaSShri Abhyankar         for (j = 0; j < bs; j++) {
44767877ebaSShri Abhyankar           for (m = 0; m < bs; m++) {
4489566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i * bs + m + shift, &row[idx]));
4499566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(bs * ajj[k] + j + shift, &col[idx]));
450a6053eceSJunchao Zhang             idx++;
45167877ebaSShri Abhyankar           }
45267877ebaSShri Abhyankar         }
45367877ebaSShri Abhyankar       }
45467877ebaSShri Abhyankar     }
455a6053eceSJunchao Zhang     mumps->irn = row;
456a6053eceSJunchao Zhang     mumps->jcn = col;
457a6053eceSJunchao Zhang     mumps->nnz = nz;
45850c845baSStefano Zampini     mumps->val = aa->a;
459127cd276SPierre Jolivet   } else if (mumps->nest_vals) PetscCall(PetscArraycpy(mumps->val, aa->a, nz)); /* MatConvertToTriples_nest_xaij() allocates mumps->val outside of MatConvertToTriples_seqbaij_seqaij(), so one needs to copy the memory */
460127cd276SPierre Jolivet   else mumps->val = aa->a;                                                      /* in the default case, mumps->val is never allocated, one just needs to update the mumps->val pointer */
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46267877ebaSShri Abhyankar }
46367877ebaSShri Abhyankar 
46466976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
465d71ae5a4SJacob Faibussowitsch {
46675480915SPierre Jolivet   const PetscInt *ai, *aj, *ajj;
467a6053eceSJunchao Zhang   PetscInt        bs;
4686497c311SBarry Smith   PetscCount      nz, rnz, i, j, k, m;
469a6053eceSJunchao Zhang   PetscMUMPSInt  *row, *col;
47075480915SPierre Jolivet   PetscScalar    *val;
47116ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa  = (Mat_SeqSBAIJ *)A->data;
47275480915SPierre Jolivet   const PetscInt  bs2 = aa->bs2, mbs = aa->mbs;
47338548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
474b94d7dedSBarry Smith   PetscBool isset, hermitian;
47538548759SBarry Smith #endif
47616ebf90aSShri Abhyankar 
47716ebf90aSShri Abhyankar   PetscFunctionBegin;
47838548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
479b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
480b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
48138548759SBarry Smith #endif
4822205254eSKarl Rupp   ai = aa->i;
4832205254eSKarl Rupp   aj = aa->j;
4849566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
48575480915SPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
4866497c311SBarry Smith     const PetscCount alloc_size = aa->nz * bs2;
487f3fa974cSJacob Faibussowitsch 
488f3fa974cSJacob Faibussowitsch     PetscCall(PetscMalloc2(alloc_size, &row, alloc_size, &col));
489a6053eceSJunchao Zhang     if (bs > 1) {
490f3fa974cSJacob Faibussowitsch       PetscCall(PetscMalloc1(alloc_size, &mumps->val_alloc));
491a6053eceSJunchao Zhang       mumps->val = mumps->val_alloc;
49275480915SPierre Jolivet     } else {
493a6053eceSJunchao Zhang       mumps->val = aa->a;
49475480915SPierre Jolivet     }
495a6053eceSJunchao Zhang     mumps->irn = row;
496a6053eceSJunchao Zhang     mumps->jcn = col;
497a6053eceSJunchao Zhang   } else {
498a6053eceSJunchao Zhang     row = mumps->irn;
499a6053eceSJunchao Zhang     col = mumps->jcn;
500a6053eceSJunchao Zhang   }
501a6053eceSJunchao Zhang   val = mumps->val;
502185f6596SHong Zhang 
50316ebf90aSShri Abhyankar   nz = 0;
504a81fe166SPierre Jolivet   if (bs > 1) {
50575480915SPierre Jolivet     for (i = 0; i < mbs; i++) {
50616ebf90aSShri Abhyankar       rnz = ai[i + 1] - ai[i];
50767877ebaSShri Abhyankar       ajj = aj + ai[i];
50875480915SPierre Jolivet       for (j = 0; j < rnz; j++) {
50975480915SPierre Jolivet         for (k = 0; k < bs; k++) {
51075480915SPierre Jolivet           for (m = 0; m < bs; m++) {
511ec4f40fdSPierre Jolivet             if (ajj[j] > i || k >= m) {
51275480915SPierre Jolivet               if (reuse == MAT_INITIAL_MATRIX) {
5139566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(i * bs + m + shift, &row[nz]));
5149566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(ajj[j] * bs + k + shift, &col[nz]));
51575480915SPierre Jolivet               }
51675480915SPierre Jolivet               val[nz++] = aa->a[(ai[i] + j) * bs2 + m + k * bs];
51775480915SPierre Jolivet             }
51875480915SPierre Jolivet           }
51975480915SPierre Jolivet         }
52075480915SPierre Jolivet       }
52175480915SPierre Jolivet     }
522a81fe166SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) {
523a81fe166SPierre Jolivet     for (i = 0; i < mbs; i++) {
524a81fe166SPierre Jolivet       rnz = ai[i + 1] - ai[i];
525a81fe166SPierre Jolivet       ajj = aj + ai[i];
526a81fe166SPierre Jolivet       for (j = 0; j < rnz; j++) {
5279566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
5289566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[nz]));
529a6053eceSJunchao Zhang         nz++;
530a81fe166SPierre Jolivet       }
531a81fe166SPierre Jolivet     }
5326497c311SBarry Smith     PetscCheck(nz == aa->nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Different numbers of nonzeros %" PetscCount_FMT " != %" PetscInt_FMT, nz, aa->nz);
533127cd276SPierre Jolivet   } else if (mumps->nest_vals)
534127cd276SPierre Jolivet     PetscCall(PetscArraycpy(mumps->val, aa->a, aa->nz)); /* bs == 1 and MAT_REUSE_MATRIX, MatConvertToTriples_nest_xaij() allocates mumps->val outside of MatConvertToTriples_seqsbaij_seqsbaij(), so one needs to copy the memory */
535127cd276SPierre Jolivet   else mumps->val = aa->a;                               /* in the default case, mumps->val is never allocated, one just needs to update the mumps->val pointer */
536a6053eceSJunchao Zhang   if (reuse == MAT_INITIAL_MATRIX) mumps->nnz = nz;
5373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53816ebf90aSShri Abhyankar }
53916ebf90aSShri Abhyankar 
54066976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
541d71ae5a4SJacob Faibussowitsch {
54267877ebaSShri Abhyankar   const PetscInt    *ai, *aj, *ajj, *adiag, M = A->rmap->n;
5436497c311SBarry Smith   PetscCount         nz, rnz, i, j;
54467877ebaSShri Abhyankar   const PetscScalar *av, *v1;
54516ebf90aSShri Abhyankar   PetscScalar       *val;
546a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
547829b1710SHong Zhang   Mat_SeqAIJ        *aa = (Mat_SeqAIJ *)A->data;
54829b521d4Sstefano_zampini   PetscBool          missing;
54938548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
550b94d7dedSBarry Smith   PetscBool hermitian, isset;
55138548759SBarry Smith #endif
55216ebf90aSShri Abhyankar 
55316ebf90aSShri Abhyankar   PetscFunctionBegin;
55438548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
555b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
556b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
55738548759SBarry Smith #endif
5589566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(A, &av));
5599371c9d4SSatish Balay   ai    = aa->i;
5609371c9d4SSatish Balay   aj    = aa->j;
56116ebf90aSShri Abhyankar   adiag = aa->diag;
5629566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal_SeqAIJ(A, &missing, NULL));
563bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
5647ee00b23SStefano Zampini     /* count nz in the upper triangular part of A */
565829b1710SHong Zhang     nz = 0;
56629b521d4Sstefano_zampini     if (missing) {
56729b521d4Sstefano_zampini       for (i = 0; i < M; i++) {
56829b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i + 1])) {
56929b521d4Sstefano_zampini           for (j = ai[i]; j < ai[i + 1]; j++) {
57029b521d4Sstefano_zampini             if (aj[j] < i) continue;
57129b521d4Sstefano_zampini             nz++;
57229b521d4Sstefano_zampini           }
57329b521d4Sstefano_zampini         } else {
57429b521d4Sstefano_zampini           nz += ai[i + 1] - adiag[i];
57529b521d4Sstefano_zampini         }
57629b521d4Sstefano_zampini       }
57729b521d4Sstefano_zampini     } else {
578829b1710SHong Zhang       for (i = 0; i < M; i++) nz += ai[i + 1] - adiag[i];
57929b521d4Sstefano_zampini     }
5809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
5819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
582a6053eceSJunchao Zhang     mumps->nnz = nz;
583a6053eceSJunchao Zhang     mumps->irn = row;
584a6053eceSJunchao Zhang     mumps->jcn = col;
585a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
586185f6596SHong Zhang 
58716ebf90aSShri Abhyankar     nz = 0;
58829b521d4Sstefano_zampini     if (missing) {
58929b521d4Sstefano_zampini       for (i = 0; i < M; i++) {
59029b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i + 1])) {
59129b521d4Sstefano_zampini           for (j = ai[i]; j < ai[i + 1]; j++) {
59229b521d4Sstefano_zampini             if (aj[j] < i) continue;
5939566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
5949566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(aj[j] + shift, &col[nz]));
59529b521d4Sstefano_zampini             val[nz] = av[j];
59629b521d4Sstefano_zampini             nz++;
59729b521d4Sstefano_zampini           }
59829b521d4Sstefano_zampini         } else {
59929b521d4Sstefano_zampini           rnz = ai[i + 1] - adiag[i];
60029b521d4Sstefano_zampini           ajj = aj + adiag[i];
60129b521d4Sstefano_zampini           v1  = av + adiag[i];
60229b521d4Sstefano_zampini           for (j = 0; j < rnz; j++) {
6039566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
6049566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[nz]));
605a6053eceSJunchao Zhang             val[nz++] = v1[j];
60629b521d4Sstefano_zampini           }
60729b521d4Sstefano_zampini         }
60829b521d4Sstefano_zampini       }
60929b521d4Sstefano_zampini     } else {
61016ebf90aSShri Abhyankar       for (i = 0; i < M; i++) {
61116ebf90aSShri Abhyankar         rnz = ai[i + 1] - adiag[i];
61267877ebaSShri Abhyankar         ajj = aj + adiag[i];
613cf3759fdSShri Abhyankar         v1  = av + adiag[i];
61467877ebaSShri Abhyankar         for (j = 0; j < rnz; j++) {
6159566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(i + shift, &row[nz]));
6169566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(ajj[j] + shift, &col[nz]));
617a6053eceSJunchao Zhang           val[nz++] = v1[j];
61816ebf90aSShri Abhyankar         }
61916ebf90aSShri Abhyankar       }
62029b521d4Sstefano_zampini     }
621397b6df1SKris Buschelman   } else {
622a6053eceSJunchao Zhang     nz  = 0;
623a6053eceSJunchao Zhang     val = mumps->val;
62429b521d4Sstefano_zampini     if (missing) {
62516ebf90aSShri Abhyankar       for (i = 0; i < M; i++) {
62629b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i + 1])) {
62729b521d4Sstefano_zampini           for (j = ai[i]; j < ai[i + 1]; j++) {
62829b521d4Sstefano_zampini             if (aj[j] < i) continue;
62929b521d4Sstefano_zampini             val[nz++] = av[j];
63029b521d4Sstefano_zampini           }
63129b521d4Sstefano_zampini         } else {
63216ebf90aSShri Abhyankar           rnz = ai[i + 1] - adiag[i];
63367877ebaSShri Abhyankar           v1  = av + adiag[i];
634ad540459SPierre Jolivet           for (j = 0; j < rnz; j++) val[nz++] = v1[j];
63516ebf90aSShri Abhyankar         }
63616ebf90aSShri Abhyankar       }
63729b521d4Sstefano_zampini     } else {
63816ebf90aSShri Abhyankar       for (i = 0; i < M; i++) {
63916ebf90aSShri Abhyankar         rnz = ai[i + 1] - adiag[i];
64016ebf90aSShri Abhyankar         v1  = av + adiag[i];
641ad540459SPierre Jolivet         for (j = 0; j < rnz; j++) val[nz++] = v1[j];
64216ebf90aSShri Abhyankar       }
64316ebf90aSShri Abhyankar     }
64429b521d4Sstefano_zampini   }
6459566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
6463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64716ebf90aSShri Abhyankar }
64816ebf90aSShri Abhyankar 
64966976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
650d71ae5a4SJacob Faibussowitsch {
651a6053eceSJunchao Zhang   const PetscInt    *ai, *aj, *bi, *bj, *garray, *ajj, *bjj;
652a6053eceSJunchao Zhang   PetscInt           bs;
6536497c311SBarry Smith   PetscCount         rstart, nz, i, j, k, m, jj, irow, countA, countB;
654a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
65516ebf90aSShri Abhyankar   const PetscScalar *av, *bv, *v1, *v2;
65616ebf90aSShri Abhyankar   PetscScalar       *val;
657397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ *)A->data;
658f4f49eeaSPierre Jolivet   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ *)mat->A->data;
659f4f49eeaSPierre Jolivet   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ *)mat->B->data;
660ec4f40fdSPierre Jolivet   const PetscInt     bs2 = aa->bs2, mbs = aa->mbs;
66138548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
662b94d7dedSBarry Smith   PetscBool hermitian, isset;
66338548759SBarry Smith #endif
66416ebf90aSShri Abhyankar 
66516ebf90aSShri Abhyankar   PetscFunctionBegin;
66638548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
667b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
668b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
66938548759SBarry Smith #endif
6709566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
67138548759SBarry Smith   rstart = A->rmap->rstart;
67238548759SBarry Smith   ai     = aa->i;
67338548759SBarry Smith   aj     = aa->j;
67438548759SBarry Smith   bi     = bb->i;
67538548759SBarry Smith   bj     = bb->j;
67638548759SBarry Smith   av     = aa->a;
67738548759SBarry Smith   bv     = bb->a;
678397b6df1SKris Buschelman 
6792205254eSKarl Rupp   garray = mat->garray;
6802205254eSKarl Rupp 
681bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
682a6053eceSJunchao Zhang     nz = (aa->nz + bb->nz) * bs2; /* just a conservative estimate */
6839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
6849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
685a6053eceSJunchao Zhang     /* can not decide the exact mumps->nnz now because of the SBAIJ */
686a6053eceSJunchao Zhang     mumps->irn = row;
687a6053eceSJunchao Zhang     mumps->jcn = col;
688a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
689397b6df1SKris Buschelman   } else {
690a6053eceSJunchao Zhang     val = mumps->val;
691397b6df1SKris Buschelman   }
692397b6df1SKris Buschelman 
6939371c9d4SSatish Balay   jj   = 0;
6949371c9d4SSatish Balay   irow = rstart;
695ec4f40fdSPierre Jolivet   for (i = 0; i < mbs; i++) {
696397b6df1SKris Buschelman     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
697397b6df1SKris Buschelman     countA = ai[i + 1] - ai[i];
698397b6df1SKris Buschelman     countB = bi[i + 1] - bi[i];
699397b6df1SKris Buschelman     bjj    = bj + bi[i];
700ec4f40fdSPierre Jolivet     v1     = av + ai[i] * bs2;
701ec4f40fdSPierre Jolivet     v2     = bv + bi[i] * bs2;
702397b6df1SKris Buschelman 
703ec4f40fdSPierre Jolivet     if (bs > 1) {
704ec4f40fdSPierre Jolivet       /* A-part */
705ec4f40fdSPierre Jolivet       for (j = 0; j < countA; j++) {
706ec4f40fdSPierre Jolivet         for (k = 0; k < bs; k++) {
707ec4f40fdSPierre Jolivet           for (m = 0; m < bs; m++) {
708ec4f40fdSPierre Jolivet             if (rstart + ajj[j] * bs > irow || k >= m) {
709ec4f40fdSPierre Jolivet               if (reuse == MAT_INITIAL_MATRIX) {
7109566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(irow + m + shift, &row[jj]));
7119566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(rstart + ajj[j] * bs + k + shift, &col[jj]));
712ec4f40fdSPierre Jolivet               }
713ec4f40fdSPierre Jolivet               val[jj++] = v1[j * bs2 + m + k * bs];
714ec4f40fdSPierre Jolivet             }
715ec4f40fdSPierre Jolivet           }
716ec4f40fdSPierre Jolivet         }
717ec4f40fdSPierre Jolivet       }
718ec4f40fdSPierre Jolivet 
719ec4f40fdSPierre Jolivet       /* B-part */
720ec4f40fdSPierre Jolivet       for (j = 0; j < countB; j++) {
721ec4f40fdSPierre Jolivet         for (k = 0; k < bs; k++) {
722ec4f40fdSPierre Jolivet           for (m = 0; m < bs; m++) {
723ec4f40fdSPierre Jolivet             if (reuse == MAT_INITIAL_MATRIX) {
7249566063dSJacob Faibussowitsch               PetscCall(PetscMUMPSIntCast(irow + m + shift, &row[jj]));
7259566063dSJacob Faibussowitsch               PetscCall(PetscMUMPSIntCast(garray[bjj[j]] * bs + k + shift, &col[jj]));
726ec4f40fdSPierre Jolivet             }
727ec4f40fdSPierre Jolivet             val[jj++] = v2[j * bs2 + m + k * bs];
728ec4f40fdSPierre Jolivet           }
729ec4f40fdSPierre Jolivet         }
730ec4f40fdSPierre Jolivet       }
731ec4f40fdSPierre Jolivet     } else {
732397b6df1SKris Buschelman       /* A-part */
733397b6df1SKris Buschelman       for (j = 0; j < countA; j++) {
734bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
7359566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
7369566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift, &col[jj]));
737397b6df1SKris Buschelman         }
73816ebf90aSShri Abhyankar         val[jj++] = v1[j];
739397b6df1SKris Buschelman       }
74016ebf90aSShri Abhyankar 
74116ebf90aSShri Abhyankar       /* B-part */
74216ebf90aSShri Abhyankar       for (j = 0; j < countB; j++) {
743bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
7449566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
7459566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift, &col[jj]));
746397b6df1SKris Buschelman         }
74716ebf90aSShri Abhyankar         val[jj++] = v2[j];
74816ebf90aSShri Abhyankar       }
74916ebf90aSShri Abhyankar     }
750ec4f40fdSPierre Jolivet     irow += bs;
751ec4f40fdSPierre Jolivet   }
7525d955bbbSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) mumps->nnz = jj;
7533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75416ebf90aSShri Abhyankar }
75516ebf90aSShri Abhyankar 
75666976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
757d71ae5a4SJacob Faibussowitsch {
75816ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj, *garray, m = A->rmap->n, *ajj, *bjj;
7596497c311SBarry Smith   PetscCount         rstart, cstart, nz, i, j, jj, irow, countA, countB;
760a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
76116ebf90aSShri Abhyankar   const PetscScalar *av, *bv, *v1, *v2;
76216ebf90aSShri Abhyankar   PetscScalar       *val;
763a3d589ffSStefano Zampini   Mat                Ad, Ao;
764a3d589ffSStefano Zampini   Mat_SeqAIJ        *aa;
765a3d589ffSStefano Zampini   Mat_SeqAIJ        *bb;
76616ebf90aSShri Abhyankar 
76716ebf90aSShri Abhyankar   PetscFunctionBegin;
7689566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A, &Ad, &Ao, &garray));
7699566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ad, &av));
7709566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ao, &bv));
771a3d589ffSStefano Zampini 
77257508eceSPierre Jolivet   aa = (Mat_SeqAIJ *)Ad->data;
77357508eceSPierre Jolivet   bb = (Mat_SeqAIJ *)Ao->data;
77438548759SBarry Smith   ai = aa->i;
77538548759SBarry Smith   aj = aa->j;
77638548759SBarry Smith   bi = bb->i;
77738548759SBarry Smith   bj = bb->j;
77816ebf90aSShri Abhyankar 
779a3d589ffSStefano Zampini   rstart = A->rmap->rstart;
7805d955bbbSStefano Zampini   cstart = A->cmap->rstart;
7812205254eSKarl Rupp 
782bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
7836497c311SBarry Smith     nz = (PetscCount)aa->nz + bb->nz; /* make sure the sum won't overflow PetscInt */
7849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
7859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
786a6053eceSJunchao Zhang     mumps->nnz = nz;
787a6053eceSJunchao Zhang     mumps->irn = row;
788a6053eceSJunchao Zhang     mumps->jcn = col;
789a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
79016ebf90aSShri Abhyankar   } else {
791a6053eceSJunchao Zhang     val = mumps->val;
79216ebf90aSShri Abhyankar   }
79316ebf90aSShri Abhyankar 
7949371c9d4SSatish Balay   jj   = 0;
7959371c9d4SSatish Balay   irow = rstart;
79616ebf90aSShri Abhyankar   for (i = 0; i < m; i++) {
79716ebf90aSShri Abhyankar     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
79816ebf90aSShri Abhyankar     countA = ai[i + 1] - ai[i];
79916ebf90aSShri Abhyankar     countB = bi[i + 1] - bi[i];
80016ebf90aSShri Abhyankar     bjj    = bj + bi[i];
80116ebf90aSShri Abhyankar     v1     = av + ai[i];
80216ebf90aSShri Abhyankar     v2     = bv + bi[i];
80316ebf90aSShri Abhyankar 
80416ebf90aSShri Abhyankar     /* A-part */
80516ebf90aSShri Abhyankar     for (j = 0; j < countA; j++) {
806bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
8079566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
8085d955bbbSStefano Zampini         PetscCall(PetscMUMPSIntCast(cstart + ajj[j] + shift, &col[jj]));
80916ebf90aSShri Abhyankar       }
81016ebf90aSShri Abhyankar       val[jj++] = v1[j];
81116ebf90aSShri Abhyankar     }
81216ebf90aSShri Abhyankar 
81316ebf90aSShri Abhyankar     /* B-part */
81416ebf90aSShri Abhyankar     for (j = 0; j < countB; j++) {
815bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
8169566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
8179566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift, &col[jj]));
81816ebf90aSShri Abhyankar       }
81916ebf90aSShri Abhyankar       val[jj++] = v2[j];
82016ebf90aSShri Abhyankar     }
82116ebf90aSShri Abhyankar     irow++;
82216ebf90aSShri Abhyankar   }
8239566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ad, &av));
8249566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ao, &bv));
8253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82616ebf90aSShri Abhyankar }
82716ebf90aSShri Abhyankar 
82866976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
829d71ae5a4SJacob Faibussowitsch {
83067877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat = (Mat_MPIBAIJ *)A->data;
831f4f49eeaSPierre Jolivet   Mat_SeqBAIJ       *aa  = (Mat_SeqBAIJ *)mat->A->data;
832f4f49eeaSPierre Jolivet   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ *)mat->B->data;
83367877ebaSShri Abhyankar   const PetscInt    *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j, *ajj, *bjj;
8345d955bbbSStefano Zampini   const PetscInt    *garray = mat->garray, mbs = mat->mbs, rstart = A->rmap->rstart, cstart = A->cmap->rstart;
83533d57670SJed Brown   const PetscInt     bs2 = mat->bs2;
836a6053eceSJunchao Zhang   PetscInt           bs;
8376497c311SBarry Smith   PetscCount         nz, i, j, k, n, jj, irow, countA, countB, idx;
838a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
83967877ebaSShri Abhyankar   const PetscScalar *av = aa->a, *bv = bb->a, *v1, *v2;
84067877ebaSShri Abhyankar   PetscScalar       *val;
84167877ebaSShri Abhyankar 
84267877ebaSShri Abhyankar   PetscFunctionBegin;
8439566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
844bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
84567877ebaSShri Abhyankar     nz = bs2 * (aa->nz + bb->nz);
8469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
8479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
848a6053eceSJunchao Zhang     mumps->nnz = nz;
849a6053eceSJunchao Zhang     mumps->irn = row;
850a6053eceSJunchao Zhang     mumps->jcn = col;
851a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
85267877ebaSShri Abhyankar   } else {
853a6053eceSJunchao Zhang     val = mumps->val;
85467877ebaSShri Abhyankar   }
85567877ebaSShri Abhyankar 
8569371c9d4SSatish Balay   jj   = 0;
8579371c9d4SSatish Balay   irow = rstart;
85867877ebaSShri Abhyankar   for (i = 0; i < mbs; i++) {
85967877ebaSShri Abhyankar     countA = ai[i + 1] - ai[i];
86067877ebaSShri Abhyankar     countB = bi[i + 1] - bi[i];
86167877ebaSShri Abhyankar     ajj    = aj + ai[i];
86267877ebaSShri Abhyankar     bjj    = bj + bi[i];
86367877ebaSShri Abhyankar     v1     = av + bs2 * ai[i];
86467877ebaSShri Abhyankar     v2     = bv + bs2 * bi[i];
86567877ebaSShri Abhyankar 
86667877ebaSShri Abhyankar     idx = 0;
86767877ebaSShri Abhyankar     /* A-part */
86867877ebaSShri Abhyankar     for (k = 0; k < countA; k++) {
86967877ebaSShri Abhyankar       for (j = 0; j < bs; j++) {
87067877ebaSShri Abhyankar         for (n = 0; n < bs; n++) {
871bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
8729566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(irow + n + shift, &row[jj]));
8735d955bbbSStefano Zampini             PetscCall(PetscMUMPSIntCast(cstart + bs * ajj[k] + j + shift, &col[jj]));
87467877ebaSShri Abhyankar           }
87567877ebaSShri Abhyankar           val[jj++] = v1[idx++];
87667877ebaSShri Abhyankar         }
87767877ebaSShri Abhyankar       }
87867877ebaSShri Abhyankar     }
87967877ebaSShri Abhyankar 
88067877ebaSShri Abhyankar     idx = 0;
88167877ebaSShri Abhyankar     /* B-part */
88267877ebaSShri Abhyankar     for (k = 0; k < countB; k++) {
88367877ebaSShri Abhyankar       for (j = 0; j < bs; j++) {
88467877ebaSShri Abhyankar         for (n = 0; n < bs; n++) {
885bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
8869566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(irow + n + shift, &row[jj]));
8879566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(bs * garray[bjj[k]] + j + shift, &col[jj]));
88867877ebaSShri Abhyankar           }
889d985c460SShri Abhyankar           val[jj++] = v2[idx++];
89067877ebaSShri Abhyankar         }
89167877ebaSShri Abhyankar       }
89267877ebaSShri Abhyankar     }
893d985c460SShri Abhyankar     irow += bs;
89467877ebaSShri Abhyankar   }
8953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89667877ebaSShri Abhyankar }
89767877ebaSShri Abhyankar 
89866976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
899d71ae5a4SJacob Faibussowitsch {
90016ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *adiag, *bi, *bj, *garray, m = A->rmap->n, *ajj, *bjj;
9016497c311SBarry Smith   PetscCount         rstart, nz, nza, nzb, i, j, jj, irow, countA, countB;
902a6053eceSJunchao Zhang   PetscMUMPSInt     *row, *col;
90316ebf90aSShri Abhyankar   const PetscScalar *av, *bv, *v1, *v2;
90416ebf90aSShri Abhyankar   PetscScalar       *val;
905a3d589ffSStefano Zampini   Mat                Ad, Ao;
906a3d589ffSStefano Zampini   Mat_SeqAIJ        *aa;
907a3d589ffSStefano Zampini   Mat_SeqAIJ        *bb;
90838548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
909b94d7dedSBarry Smith   PetscBool hermitian, isset;
91038548759SBarry Smith #endif
91116ebf90aSShri Abhyankar 
91216ebf90aSShri Abhyankar   PetscFunctionBegin;
91338548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
914b94d7dedSBarry Smith   PetscCall(MatIsHermitianKnown(A, &isset, &hermitian));
915b94d7dedSBarry Smith   PetscCheck(!isset || !hermitian, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MUMPS does not support Hermitian symmetric matrices for Choleksy");
91638548759SBarry Smith #endif
9179566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A, &Ad, &Ao, &garray));
9189566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ad, &av));
9199566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ao, &bv));
920a3d589ffSStefano Zampini 
92157508eceSPierre Jolivet   aa    = (Mat_SeqAIJ *)Ad->data;
92257508eceSPierre Jolivet   bb    = (Mat_SeqAIJ *)Ao->data;
92338548759SBarry Smith   ai    = aa->i;
92438548759SBarry Smith   aj    = aa->j;
92538548759SBarry Smith   adiag = aa->diag;
92638548759SBarry Smith   bi    = bb->i;
92738548759SBarry Smith   bj    = bb->j;
9282205254eSKarl Rupp 
92916ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
93016ebf90aSShri Abhyankar 
931bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
932e0bace9bSHong Zhang     nza = 0; /* num of upper triangular entries in mat->A, including diagonals */
933e0bace9bSHong Zhang     nzb = 0; /* num of upper triangular entries in mat->B */
93416ebf90aSShri Abhyankar     for (i = 0; i < m; i++) {
935e0bace9bSHong Zhang       nza += (ai[i + 1] - adiag[i]);
93616ebf90aSShri Abhyankar       countB = bi[i + 1] - bi[i];
93716ebf90aSShri Abhyankar       bjj    = bj + bi[i];
938e0bace9bSHong Zhang       for (j = 0; j < countB; j++) {
939e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
940e0bace9bSHong Zhang       }
941e0bace9bSHong Zhang     }
94216ebf90aSShri Abhyankar 
943e0bace9bSHong Zhang     nz = nza + nzb; /* total nz of upper triangular part of mat */
9449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &row, nz, &col));
9459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &val));
946a6053eceSJunchao Zhang     mumps->nnz = nz;
947a6053eceSJunchao Zhang     mumps->irn = row;
948a6053eceSJunchao Zhang     mumps->jcn = col;
949a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
95016ebf90aSShri Abhyankar   } else {
951a6053eceSJunchao Zhang     val = mumps->val;
95216ebf90aSShri Abhyankar   }
95316ebf90aSShri Abhyankar 
9549371c9d4SSatish Balay   jj   = 0;
9559371c9d4SSatish Balay   irow = rstart;
95616ebf90aSShri Abhyankar   for (i = 0; i < m; i++) {
95716ebf90aSShri Abhyankar     ajj    = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */
95816ebf90aSShri Abhyankar     v1     = av + adiag[i];
95916ebf90aSShri Abhyankar     countA = ai[i + 1] - adiag[i];
96016ebf90aSShri Abhyankar     countB = bi[i + 1] - bi[i];
96116ebf90aSShri Abhyankar     bjj    = bj + bi[i];
96216ebf90aSShri Abhyankar     v2     = bv + bi[i];
96316ebf90aSShri Abhyankar 
96416ebf90aSShri Abhyankar     /* A-part */
96516ebf90aSShri Abhyankar     for (j = 0; j < countA; j++) {
966bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
9679566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
9689566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift, &col[jj]));
96916ebf90aSShri Abhyankar       }
97016ebf90aSShri Abhyankar       val[jj++] = v1[j];
97116ebf90aSShri Abhyankar     }
97216ebf90aSShri Abhyankar 
97316ebf90aSShri Abhyankar     /* B-part */
97416ebf90aSShri Abhyankar     for (j = 0; j < countB; j++) {
97516ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
976bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
9779566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift, &row[jj]));
9789566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift, &col[jj]));
97916ebf90aSShri Abhyankar         }
98016ebf90aSShri Abhyankar         val[jj++] = v2[j];
98116ebf90aSShri Abhyankar       }
982397b6df1SKris Buschelman     }
983397b6df1SKris Buschelman     irow++;
984397b6df1SKris Buschelman   }
9859566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ad, &av));
9869566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ao, &bv));
9873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
988397b6df1SKris Buschelman }
989397b6df1SKris Buschelman 
990d2a308c1SPierre Jolivet static PetscErrorCode MatConvertToTriples_diagonal_xaij(Mat A, PETSC_UNUSED PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
991c3e1b152SPierre Jolivet {
992c3e1b152SPierre Jolivet   const PetscScalar *av;
993c3e1b152SPierre Jolivet   const PetscInt     M = A->rmap->n;
9946497c311SBarry Smith   PetscCount         i;
995c3e1b152SPierre Jolivet   PetscMUMPSInt     *row, *col;
996c3e1b152SPierre Jolivet   Vec                v;
997c3e1b152SPierre Jolivet 
998c3e1b152SPierre Jolivet   PetscFunctionBegin;
999c3e1b152SPierre Jolivet   PetscCall(MatDiagonalGetDiagonal(A, &v));
1000c3e1b152SPierre Jolivet   PetscCall(VecGetArrayRead(v, &av));
1001c3e1b152SPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
1002c3e1b152SPierre Jolivet     PetscCall(PetscMalloc2(M, &row, M, &col));
1003c3e1b152SPierre Jolivet     for (i = 0; i < M; i++) {
1004c3e1b152SPierre Jolivet       PetscCall(PetscMUMPSIntCast(i + A->rmap->rstart, &row[i]));
1005c3e1b152SPierre Jolivet       col[i] = row[i];
1006c3e1b152SPierre Jolivet     }
1007c3e1b152SPierre Jolivet     mumps->val = (PetscScalar *)av;
1008c3e1b152SPierre Jolivet     mumps->irn = row;
1009c3e1b152SPierre Jolivet     mumps->jcn = col;
1010c3e1b152SPierre Jolivet     mumps->nnz = M;
1011127cd276SPierre Jolivet   } else if (mumps->nest_vals) PetscCall(PetscArraycpy(mumps->val, av, M)); /* MatConvertToTriples_nest_xaij() allocates mumps->val outside of MatConvertToTriples_diagonal_xaij(), so one needs to copy the memory */
1012127cd276SPierre Jolivet   else mumps->val = (PetscScalar *)av;                                      /* in the default case, mumps->val is never allocated, one just needs to update the mumps->val pointer */
1013c3e1b152SPierre Jolivet   PetscCall(VecRestoreArrayRead(v, &av));
1014c3e1b152SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
1015c3e1b152SPierre Jolivet }
1016c3e1b152SPierre Jolivet 
1017d2a308c1SPierre Jolivet static PetscErrorCode MatConvertToTriples_dense_xaij(Mat A, PETSC_UNUSED PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
10184b9405b2SPierre Jolivet {
10194b9405b2SPierre Jolivet   PetscScalar   *v;
10204b9405b2SPierre Jolivet   const PetscInt m = A->rmap->n, N = A->cmap->N;
10214b9405b2SPierre Jolivet   PetscInt       lda;
10226497c311SBarry Smith   PetscCount     i, j;
10234b9405b2SPierre Jolivet   PetscMUMPSInt *row, *col;
10244b9405b2SPierre Jolivet 
10254b9405b2SPierre Jolivet   PetscFunctionBegin;
10264b9405b2SPierre Jolivet   PetscCall(MatDenseGetArray(A, &v));
10274b9405b2SPierre Jolivet   PetscCall(MatDenseGetLDA(A, &lda));
10284b9405b2SPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
10294b9405b2SPierre Jolivet     PetscCall(PetscMalloc2(m * N, &row, m * N, &col));
10304b9405b2SPierre Jolivet     for (i = 0; i < m; i++) {
10314b9405b2SPierre Jolivet       col[i] = 0;
10324b9405b2SPierre Jolivet       PetscCall(PetscMUMPSIntCast(i + A->rmap->rstart, &row[i]));
10334b9405b2SPierre Jolivet     }
10344b9405b2SPierre Jolivet     for (j = 1; j < N; j++) {
10354b9405b2SPierre Jolivet       for (i = 0; i < m; i++) PetscCall(PetscMUMPSIntCast(j, col + i + m * j));
10364b9405b2SPierre Jolivet       PetscCall(PetscArraycpy(row + m * j, row + m * (j - 1), m));
10374b9405b2SPierre Jolivet     }
10384b9405b2SPierre Jolivet     if (lda == m) mumps->val = v;
10394b9405b2SPierre Jolivet     else {
10404b9405b2SPierre Jolivet       PetscCall(PetscMalloc1(m * N, &mumps->val));
10414b9405b2SPierre Jolivet       mumps->val_alloc = mumps->val;
10424b9405b2SPierre Jolivet       for (j = 0; j < N; j++) PetscCall(PetscArraycpy(mumps->val + m * j, v + lda * j, m));
10434b9405b2SPierre Jolivet     }
10444b9405b2SPierre Jolivet     mumps->irn = row;
10454b9405b2SPierre Jolivet     mumps->jcn = col;
10464b9405b2SPierre Jolivet     mumps->nnz = m * N;
10474b9405b2SPierre Jolivet   } else {
10484b9405b2SPierre Jolivet     if (lda == m && !mumps->nest_vals) mumps->val = v;
10494b9405b2SPierre Jolivet     else {
10504b9405b2SPierre Jolivet       for (j = 0; j < N; j++) PetscCall(PetscArraycpy(mumps->val + m * j, v + lda * j, m));
10514b9405b2SPierre Jolivet     }
10524b9405b2SPierre Jolivet   }
10534b9405b2SPierre Jolivet   PetscCall(MatDenseRestoreArray(A, &v));
10544b9405b2SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10554b9405b2SPierre Jolivet }
10564b9405b2SPierre Jolivet 
105766976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertToTriples_nest_xaij(Mat A, PetscInt shift, MatReuse reuse, Mat_MUMPS *mumps)
10589d0448ceSStefano Zampini {
10599d0448ceSStefano Zampini   Mat     **mats;
10609d0448ceSStefano Zampini   PetscInt  nr, nc;
10619d0448ceSStefano Zampini   PetscBool chol = mumps->sym ? PETSC_TRUE : PETSC_FALSE;
10629d0448ceSStefano Zampini 
10639d0448ceSStefano Zampini   PetscFunctionBegin;
10649d0448ceSStefano Zampini   PetscCall(MatNestGetSubMats(A, &nr, &nc, &mats));
10659d0448ceSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
10669d0448ceSStefano Zampini     PetscMUMPSInt *irns, *jcns;
10679d0448ceSStefano Zampini     PetscScalar   *vals;
10686497c311SBarry Smith     PetscCount     totnnz, cumnnz, maxnnz;
10699d0448ceSStefano Zampini     PetscInt      *pjcns_w;
10709d0448ceSStefano Zampini     IS            *rows, *cols;
10719d0448ceSStefano Zampini     PetscInt     **rows_idx, **cols_idx;
10729d0448ceSStefano Zampini 
10739d0448ceSStefano Zampini     cumnnz = 0;
10749d0448ceSStefano Zampini     maxnnz = 0;
10755d955bbbSStefano Zampini     PetscCall(PetscMalloc2(nr * nc + 1, &mumps->nest_vals_start, nr * nc, &mumps->nest_convert_to_triples));
10769d0448ceSStefano Zampini     for (PetscInt r = 0; r < nr; r++) {
10779d0448ceSStefano Zampini       for (PetscInt c = 0; c < nc; c++) {
10789d0448ceSStefano Zampini         Mat sub = mats[r][c];
10799d0448ceSStefano Zampini 
10809d0448ceSStefano Zampini         mumps->nest_convert_to_triples[r * nc + c] = NULL;
10819d0448ceSStefano Zampini         if (chol && c < r) continue; /* skip lower-triangular block for Cholesky */
10829d0448ceSStefano Zampini         if (sub) {
10839d0448ceSStefano Zampini           PetscErrorCode (*convert_to_triples)(Mat, PetscInt, MatReuse, Mat_MUMPS *) = NULL;
10844b9405b2SPierre Jolivet           PetscBool isSeqAIJ, isMPIAIJ, isSeqBAIJ, isMPIBAIJ, isSeqSBAIJ, isMPISBAIJ, isTrans, isHTrans = PETSC_FALSE, isDiag, isDense;
10859d0448ceSStefano Zampini           MatInfo   info;
10869d0448ceSStefano Zampini 
10875d955bbbSStefano Zampini           PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATTRANSPOSEVIRTUAL, &isTrans));
10885d955bbbSStefano Zampini           if (isTrans) PetscCall(MatTransposeGetMat(sub, &sub));
10895d955bbbSStefano Zampini           else {
10905d955bbbSStefano Zampini             PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATHERMITIANTRANSPOSEVIRTUAL, &isHTrans));
10915d955bbbSStefano Zampini             if (isHTrans) PetscCall(MatHermitianTransposeGetMat(sub, &sub));
10925d955bbbSStefano Zampini           }
10939d0448ceSStefano Zampini           PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATSEQAIJ, &isSeqAIJ));
10949d0448ceSStefano Zampini           PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATMPIAIJ, &isMPIAIJ));
10959d0448ceSStefano Zampini           PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATSEQBAIJ, &isSeqBAIJ));
10969d0448ceSStefano Zampini           PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATMPIBAIJ, &isMPIBAIJ));
10979d0448ceSStefano Zampini           PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATSEQSBAIJ, &isSeqSBAIJ));
10989d0448ceSStefano Zampini           PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATMPISBAIJ, &isMPISBAIJ));
1099c3e1b152SPierre Jolivet           PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATDIAGONAL, &isDiag));
11004b9405b2SPierre Jolivet           PetscCall(PetscObjectTypeCompareAny((PetscObject)sub, &isDense, MATSEQDENSE, MATMPIDENSE, NULL));
11019d0448ceSStefano Zampini 
11029d0448ceSStefano Zampini           if (chol) {
11039d0448ceSStefano Zampini             if (r == c) {
11049d0448ceSStefano Zampini               if (isSeqAIJ) convert_to_triples = MatConvertToTriples_seqaij_seqsbaij;
11059d0448ceSStefano Zampini               else if (isMPIAIJ) convert_to_triples = MatConvertToTriples_mpiaij_mpisbaij;
11069d0448ceSStefano Zampini               else if (isSeqSBAIJ) convert_to_triples = MatConvertToTriples_seqsbaij_seqsbaij;
11079d0448ceSStefano Zampini               else if (isMPISBAIJ) convert_to_triples = MatConvertToTriples_mpisbaij_mpisbaij;
1108c3e1b152SPierre Jolivet               else if (isDiag) convert_to_triples = MatConvertToTriples_diagonal_xaij;
11094b9405b2SPierre Jolivet               else if (isDense) convert_to_triples = MatConvertToTriples_dense_xaij;
11109d0448ceSStefano Zampini             } else {
11119d0448ceSStefano Zampini               if (isSeqAIJ) convert_to_triples = MatConvertToTriples_seqaij_seqaij;
11129d0448ceSStefano Zampini               else if (isMPIAIJ) convert_to_triples = MatConvertToTriples_mpiaij_mpiaij;
11139d0448ceSStefano Zampini               else if (isSeqBAIJ) convert_to_triples = MatConvertToTriples_seqbaij_seqaij;
11149d0448ceSStefano Zampini               else if (isMPIBAIJ) convert_to_triples = MatConvertToTriples_mpibaij_mpiaij;
1115c3e1b152SPierre Jolivet               else if (isDiag) convert_to_triples = MatConvertToTriples_diagonal_xaij;
11164b9405b2SPierre Jolivet               else if (isDense) convert_to_triples = MatConvertToTriples_dense_xaij;
11179d0448ceSStefano Zampini             }
11189d0448ceSStefano Zampini           } else {
11199d0448ceSStefano Zampini             if (isSeqAIJ) convert_to_triples = MatConvertToTriples_seqaij_seqaij;
11209d0448ceSStefano Zampini             else if (isMPIAIJ) convert_to_triples = MatConvertToTriples_mpiaij_mpiaij;
11219d0448ceSStefano Zampini             else if (isSeqBAIJ) convert_to_triples = MatConvertToTriples_seqbaij_seqaij;
11229d0448ceSStefano Zampini             else if (isMPIBAIJ) convert_to_triples = MatConvertToTriples_mpibaij_mpiaij;
1123c3e1b152SPierre Jolivet             else if (isDiag) convert_to_triples = MatConvertToTriples_diagonal_xaij;
11244b9405b2SPierre Jolivet             else if (isDense) convert_to_triples = MatConvertToTriples_dense_xaij;
11259d0448ceSStefano Zampini           }
11269d0448ceSStefano Zampini           PetscCheck(convert_to_triples, PetscObjectComm((PetscObject)sub), PETSC_ERR_SUP, "Not for block of type %s", ((PetscObject)sub)->type_name);
11279d0448ceSStefano Zampini           mumps->nest_convert_to_triples[r * nc + c] = convert_to_triples;
11289d0448ceSStefano Zampini           PetscCall(MatGetInfo(sub, MAT_LOCAL, &info));
11296497c311SBarry Smith           cumnnz += (PetscCount)info.nz_used; /* can be overestimated for Cholesky */
11309d0448ceSStefano Zampini           maxnnz = PetscMax(maxnnz, info.nz_used);
11319d0448ceSStefano Zampini         }
11329d0448ceSStefano Zampini       }
11339d0448ceSStefano Zampini     }
11349d0448ceSStefano Zampini 
11359d0448ceSStefano Zampini     /* Allocate total COO */
11369d0448ceSStefano Zampini     totnnz = cumnnz;
11379d0448ceSStefano Zampini     PetscCall(PetscMalloc2(totnnz, &irns, totnnz, &jcns));
11389d0448ceSStefano Zampini     PetscCall(PetscMalloc1(totnnz, &vals));
11399d0448ceSStefano Zampini 
11409d0448ceSStefano Zampini     /* Handle rows and column maps
11419d0448ceSStefano Zampini        We directly map rows and use an SF for the columns */
11429d0448ceSStefano Zampini     PetscCall(PetscMalloc4(nr, &rows, nc, &cols, nr, &rows_idx, nc, &cols_idx));
11439d0448ceSStefano Zampini     PetscCall(MatNestGetISs(A, rows, cols));
11449d0448ceSStefano Zampini     for (PetscInt r = 0; r < nr; r++) PetscCall(ISGetIndices(rows[r], (const PetscInt **)&rows_idx[r]));
11459d0448ceSStefano Zampini     for (PetscInt c = 0; c < nc; c++) PetscCall(ISGetIndices(cols[c], (const PetscInt **)&cols_idx[c]));
11469d0448ceSStefano Zampini     if (PetscDefined(USE_64BIT_INDICES)) PetscCall(PetscMalloc1(maxnnz, &pjcns_w));
11475d955bbbSStefano Zampini     else (void)maxnnz;
11489d0448ceSStefano Zampini 
11499d0448ceSStefano Zampini     cumnnz = 0;
11509d0448ceSStefano Zampini     for (PetscInt r = 0; r < nr; r++) {
11519d0448ceSStefano Zampini       for (PetscInt c = 0; c < nc; c++) {
11529d0448ceSStefano Zampini         Mat             sub  = mats[r][c];
11539d0448ceSStefano Zampini         const PetscInt *ridx = rows_idx[r];
11545d955bbbSStefano Zampini         const PetscInt *cidx = cols_idx[c];
11559d0448ceSStefano Zampini         PetscInt        rst;
11569d0448ceSStefano Zampini         PetscSF         csf;
11575d955bbbSStefano Zampini         PetscBool       isTrans, isHTrans = PETSC_FALSE, swap;
11585d955bbbSStefano Zampini         PetscLayout     cmap;
11596497c311SBarry Smith         PetscInt        innz;
11609d0448ceSStefano Zampini 
11619d0448ceSStefano Zampini         mumps->nest_vals_start[r * nc + c] = cumnnz;
11629d0448ceSStefano Zampini         if (!mumps->nest_convert_to_triples[r * nc + c]) continue;
11639d0448ceSStefano Zampini 
11645d955bbbSStefano Zampini         /* Extract inner blocks if needed */
11655d955bbbSStefano Zampini         PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATTRANSPOSEVIRTUAL, &isTrans));
11665d955bbbSStefano Zampini         if (isTrans) PetscCall(MatTransposeGetMat(sub, &sub));
11675d955bbbSStefano Zampini         else {
11685d955bbbSStefano Zampini           PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATHERMITIANTRANSPOSEVIRTUAL, &isHTrans));
11695d955bbbSStefano Zampini           if (isHTrans) PetscCall(MatHermitianTransposeGetMat(sub, &sub));
11705d955bbbSStefano Zampini         }
11715d955bbbSStefano Zampini         swap = (PetscBool)(isTrans || isHTrans);
11725d955bbbSStefano Zampini 
11735d955bbbSStefano Zampini         /* Get column layout to map off-process columns */
11745d955bbbSStefano Zampini         PetscCall(MatGetLayouts(sub, NULL, &cmap));
11755d955bbbSStefano Zampini 
11765d955bbbSStefano Zampini         /* Get row start to map on-process rows */
11775d955bbbSStefano Zampini         PetscCall(MatGetOwnershipRange(sub, &rst, NULL));
11785d955bbbSStefano Zampini 
11799d0448ceSStefano Zampini         /* Directly use the mumps datastructure and use C ordering for now */
11809d0448ceSStefano Zampini         PetscCall((*mumps->nest_convert_to_triples[r * nc + c])(sub, 0, MAT_INITIAL_MATRIX, mumps));
11819d0448ceSStefano Zampini 
11825d955bbbSStefano Zampini         /* Swap the role of rows and columns indices for transposed blocks
11835d955bbbSStefano Zampini            since we need values with global final ordering */
11845d955bbbSStefano Zampini         if (swap) {
11855d955bbbSStefano Zampini           cidx = rows_idx[r];
11865d955bbbSStefano Zampini           ridx = cols_idx[c];
11879d0448ceSStefano Zampini         }
11889d0448ceSStefano Zampini 
11895d955bbbSStefano Zampini         /* Communicate column indices
11905d955bbbSStefano Zampini            This could have been done with a single SF but it would have complicated the code a lot.
11915d955bbbSStefano Zampini            But since we do it only once, we pay the price of setting up an SF for each block */
11925d955bbbSStefano Zampini         if (PetscDefined(USE_64BIT_INDICES)) {
11935d955bbbSStefano Zampini           for (PetscInt k = 0; k < mumps->nnz; k++) pjcns_w[k] = mumps->jcn[k];
1194f4f49eeaSPierre Jolivet         } else pjcns_w = (PetscInt *)mumps->jcn; /* This cast is needed only to silence warnings for 64bit integers builds */
11959d0448ceSStefano Zampini         PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &csf));
11966497c311SBarry Smith         PetscCall(PetscIntCast(mumps->nnz, &innz));
11976497c311SBarry Smith         PetscCall(PetscSFSetGraphLayout(csf, cmap, innz, NULL, PETSC_OWN_POINTER, pjcns_w));
11985d955bbbSStefano Zampini         PetscCall(PetscSFBcastBegin(csf, MPIU_INT, cidx, pjcns_w, MPI_REPLACE));
11995d955bbbSStefano Zampini         PetscCall(PetscSFBcastEnd(csf, MPIU_INT, cidx, pjcns_w, MPI_REPLACE));
12009d0448ceSStefano Zampini         PetscCall(PetscSFDestroy(&csf));
12019d0448ceSStefano Zampini 
12025d955bbbSStefano Zampini         /* Import indices: use direct map for rows and mapped indices for columns */
12035d955bbbSStefano Zampini         if (swap) {
12045d955bbbSStefano Zampini           for (PetscInt k = 0; k < mumps->nnz; k++) {
12055d955bbbSStefano Zampini             PetscCall(PetscMUMPSIntCast(ridx[mumps->irn[k] - rst] + shift, &jcns[cumnnz + k]));
12065d955bbbSStefano Zampini             PetscCall(PetscMUMPSIntCast(pjcns_w[k] + shift, &irns[cumnnz + k]));
12075d955bbbSStefano Zampini           }
12085d955bbbSStefano Zampini         } else {
12095d955bbbSStefano Zampini           for (PetscInt k = 0; k < mumps->nnz; k++) {
12105d955bbbSStefano Zampini             PetscCall(PetscMUMPSIntCast(ridx[mumps->irn[k] - rst] + shift, &irns[cumnnz + k]));
12115d955bbbSStefano Zampini             PetscCall(PetscMUMPSIntCast(pjcns_w[k] + shift, &jcns[cumnnz + k]));
12125d955bbbSStefano Zampini           }
12135d955bbbSStefano Zampini         }
12145d955bbbSStefano Zampini 
12155d955bbbSStefano Zampini         /* Import values to full COO */
12165d955bbbSStefano Zampini         PetscCall(PetscArraycpy(vals + cumnnz, mumps->val, mumps->nnz));
121750c845baSStefano Zampini         if (isHTrans) { /* conjugate the entries */
121850c845baSStefano Zampini           PetscScalar *v = vals + cumnnz;
121950c845baSStefano Zampini           for (PetscInt k = 0; k < mumps->nnz; k++) v[k] = PetscConj(v[k]);
122050c845baSStefano Zampini         }
12219d0448ceSStefano Zampini 
12229d0448ceSStefano Zampini         /* Shift new starting point and sanity check */
12239d0448ceSStefano Zampini         cumnnz += mumps->nnz;
12246497c311SBarry Smith         PetscCheck(cumnnz <= totnnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected number of nonzeros %" PetscCount_FMT " != %" PetscCount_FMT, cumnnz, totnnz);
12259d0448ceSStefano Zampini 
12269d0448ceSStefano Zampini         /* Free scratch memory */
12279d0448ceSStefano Zampini         PetscCall(PetscFree2(mumps->irn, mumps->jcn));
12289d0448ceSStefano Zampini         PetscCall(PetscFree(mumps->val_alloc));
12299d0448ceSStefano Zampini         mumps->val = NULL;
12309d0448ceSStefano Zampini         mumps->nnz = 0;
12319d0448ceSStefano Zampini       }
12329d0448ceSStefano Zampini     }
12339d0448ceSStefano Zampini     if (PetscDefined(USE_64BIT_INDICES)) PetscCall(PetscFree(pjcns_w));
12349d0448ceSStefano Zampini     for (PetscInt r = 0; r < nr; r++) PetscCall(ISRestoreIndices(rows[r], (const PetscInt **)&rows_idx[r]));
12359d0448ceSStefano Zampini     for (PetscInt c = 0; c < nc; c++) PetscCall(ISRestoreIndices(cols[c], (const PetscInt **)&cols_idx[c]));
12369d0448ceSStefano Zampini     PetscCall(PetscFree4(rows, cols, rows_idx, cols_idx));
12376497c311SBarry Smith     if (!chol) PetscCheck(cumnnz == totnnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Different number of nonzeros %" PetscCount_FMT " != %" PetscCount_FMT, cumnnz, totnnz);
12385d955bbbSStefano Zampini     mumps->nest_vals_start[nr * nc] = cumnnz;
12399d0448ceSStefano Zampini 
12409d0448ceSStefano Zampini     /* Set pointers for final MUMPS data structure */
12419d0448ceSStefano Zampini     mumps->nest_vals = vals;
12429d0448ceSStefano Zampini     mumps->val_alloc = NULL; /* do not use val_alloc since it may be reallocated with the OMP callpath */
12439d0448ceSStefano Zampini     mumps->val       = vals;
12449d0448ceSStefano Zampini     mumps->irn       = irns;
12459d0448ceSStefano Zampini     mumps->jcn       = jcns;
12469d0448ceSStefano Zampini     mumps->nnz       = cumnnz;
12479d0448ceSStefano Zampini   } else {
12489d0448ceSStefano Zampini     PetscScalar *oval = mumps->nest_vals;
12499d0448ceSStefano Zampini     for (PetscInt r = 0; r < nr; r++) {
12509d0448ceSStefano Zampini       for (PetscInt c = 0; c < nc; c++) {
12515d955bbbSStefano Zampini         PetscBool isTrans, isHTrans = PETSC_FALSE;
12525d955bbbSStefano Zampini         Mat       sub  = mats[r][c];
12535d955bbbSStefano Zampini         PetscInt  midx = r * nc + c;
12545d955bbbSStefano Zampini 
12555d955bbbSStefano Zampini         if (!mumps->nest_convert_to_triples[midx]) continue;
12565d955bbbSStefano Zampini         PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATTRANSPOSEVIRTUAL, &isTrans));
12575d955bbbSStefano Zampini         if (isTrans) PetscCall(MatTransposeGetMat(sub, &sub));
12585d955bbbSStefano Zampini         else {
12595d955bbbSStefano Zampini           PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATHERMITIANTRANSPOSEVIRTUAL, &isHTrans));
12605d955bbbSStefano Zampini           if (isHTrans) PetscCall(MatHermitianTransposeGetMat(sub, &sub));
12615d955bbbSStefano Zampini         }
12625d955bbbSStefano Zampini         mumps->val = oval + mumps->nest_vals_start[midx];
12635d955bbbSStefano Zampini         PetscCall((*mumps->nest_convert_to_triples[midx])(sub, shift, MAT_REUSE_MATRIX, mumps));
12645d955bbbSStefano Zampini         if (isHTrans) {
12656497c311SBarry Smith           PetscCount nnz = mumps->nest_vals_start[midx + 1] - mumps->nest_vals_start[midx];
12666497c311SBarry Smith           for (PetscCount k = 0; k < nnz; k++) mumps->val[k] = PetscConj(mumps->val[k]);
12675d955bbbSStefano Zampini         }
12689d0448ceSStefano Zampini       }
12699d0448ceSStefano Zampini     }
12709d0448ceSStefano Zampini     mumps->val = oval;
12719d0448ceSStefano Zampini   }
12729d0448ceSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
12739d0448ceSStefano Zampini }
12749d0448ceSStefano Zampini 
127566976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_MUMPS(Mat A)
1276d71ae5a4SJacob Faibussowitsch {
1277a6053eceSJunchao Zhang   Mat_MUMPS *mumps = (Mat_MUMPS *)A->data;
1278b24902e0SBarry Smith 
1279397b6df1SKris Buschelman   PetscFunctionBegin;
12809566063dSJacob Faibussowitsch   PetscCall(PetscFree2(mumps->id.sol_loc, mumps->id.isol_loc));
12819566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mumps->scat_rhs));
12829566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mumps->scat_sol));
12839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mumps->b_seq));
12849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mumps->x_seq));
12859566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.perm_in));
12869566063dSJacob Faibussowitsch   PetscCall(PetscFree2(mumps->irn, mumps->jcn));
12879566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->val_alloc));
12889566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->info));
1289413bcc21SPierre Jolivet   PetscCall(PetscFree(mumps->ICNTL_pre));
1290413bcc21SPierre Jolivet   PetscCall(PetscFree(mumps->CNTL_pre));
12919566063dSJacob Faibussowitsch   PetscCall(MatMumpsResetSchur_Private(mumps));
1292413bcc21SPierre Jolivet   if (mumps->id.job != JOB_NULL) { /* cannot call PetscMUMPS_c() if JOB_INIT has never been called for this instance */
1293a5e57a09SHong Zhang     mumps->id.job = JOB_END;
12943ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
12959261f6e4SBarry Smith     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in termination: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
1296413bcc21SPierre Jolivet     if (mumps->mumps_comm != MPI_COMM_NULL) {
1297413bcc21SPierre Jolivet       if (PetscDefined(HAVE_OPENMP_SUPPORT) && mumps->use_petsc_omp_support) PetscCallMPI(MPI_Comm_free(&mumps->mumps_comm));
1298413bcc21SPierre Jolivet       else PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)A), &mumps->mumps_comm));
1299413bcc21SPierre Jolivet     }
1300413bcc21SPierre Jolivet   }
13013ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
130267602552SJunchao Zhang   if (mumps->use_petsc_omp_support) {
13039566063dSJacob Faibussowitsch     PetscCall(PetscOmpCtrlDestroy(&mumps->omp_ctrl));
13049566063dSJacob Faibussowitsch     PetscCall(PetscFree2(mumps->rhs_loc, mumps->rhs_recvbuf));
13059566063dSJacob Faibussowitsch     PetscCall(PetscFree3(mumps->rhs_nrow, mumps->rhs_recvcounts, mumps->rhs_disps));
130667602552SJunchao Zhang   }
13073ab56b82SJunchao Zhang #endif
13089566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->ia_alloc));
13099566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->ja_alloc));
13109566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->recvcount));
13119566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->reqs));
13129566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->irhs_loc));
13139d0448ceSStefano Zampini   PetscCall(PetscFree2(mumps->nest_vals_start, mumps->nest_convert_to_triples));
13149d0448ceSStefano Zampini   PetscCall(PetscFree(mumps->nest_vals));
13159566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1316bf0cc555SLisandro Dalcin 
131797969023SHong Zhang   /* clear composed functions */
13189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
13199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorSetSchurIS_C", NULL));
13209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorCreateSchurComplement_C", NULL));
13219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsSetIcntl_C", NULL));
13229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetIcntl_C", NULL));
13239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsSetCntl_C", NULL));
13249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetCntl_C", NULL));
13259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInfo_C", NULL));
13269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInfog_C", NULL));
13279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetRinfo_C", NULL));
13289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetRinfog_C", NULL));
13295c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetNullPivots_C", NULL));
13309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInverse_C", NULL));
13319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMumpsGetInverseTranspose_C", NULL));
13323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1333397b6df1SKris Buschelman }
1334397b6df1SKris Buschelman 
133567602552SJunchao 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. */
1336d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMumpsSetUpDistRHSInfo(Mat A, PetscInt nrhs, const PetscScalar *array)
1337d71ae5a4SJacob Faibussowitsch {
133867602552SJunchao Zhang   Mat_MUMPS        *mumps   = (Mat_MUMPS *)A->data;
133967602552SJunchao Zhang   const PetscMPIInt ompsize = mumps->omp_comm_size;
134067602552SJunchao Zhang   PetscInt          i, m, M, rstart;
134167602552SJunchao Zhang 
134267602552SJunchao Zhang   PetscFunctionBegin;
13439566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, NULL));
13449566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
134508401ef6SPierre Jolivet   PetscCheck(M <= PETSC_MUMPS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "PetscInt too long for PetscMUMPSInt");
134667602552SJunchao Zhang   if (ompsize == 1) {
134767602552SJunchao Zhang     if (!mumps->irhs_loc) {
13486497c311SBarry Smith       mumps->nloc_rhs = (PetscMUMPSInt)m;
13499566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &mumps->irhs_loc));
13509566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
13516497c311SBarry Smith       for (i = 0; i < m; i++) PetscCall(PetscMUMPSIntCast(rstart + i + 1, &mumps->irhs_loc[i])); /* use 1-based indices */
135267602552SJunchao Zhang     }
135367602552SJunchao Zhang     mumps->id.rhs_loc = (MumpsScalar *)array;
135467602552SJunchao Zhang   } else {
135567602552SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
135667602552SJunchao Zhang     const PetscInt *ranges;
135767602552SJunchao Zhang     PetscMPIInt     j, k, sendcount, *petsc_ranks, *omp_ranks;
135867602552SJunchao Zhang     MPI_Group       petsc_group, omp_group;
135967602552SJunchao Zhang     PetscScalar    *recvbuf = NULL;
136067602552SJunchao Zhang 
136167602552SJunchao Zhang     if (mumps->is_omp_master) {
136267602552SJunchao Zhang       /* Lazily initialize the omp stuff for distributed rhs */
136367602552SJunchao Zhang       if (!mumps->irhs_loc) {
13649566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(ompsize, &omp_ranks, ompsize, &petsc_ranks));
13659566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(ompsize, &mumps->rhs_nrow, ompsize, &mumps->rhs_recvcounts, ompsize, &mumps->rhs_disps));
13669566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_group(mumps->petsc_comm, &petsc_group));
13679566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_group(mumps->omp_comm, &omp_group));
136867602552SJunchao Zhang         for (j = 0; j < ompsize; j++) omp_ranks[j] = j;
13699566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_translate_ranks(omp_group, ompsize, omp_ranks, petsc_group, petsc_ranks));
137067602552SJunchao Zhang 
137167602552SJunchao Zhang         /* Populate mumps->irhs_loc[], rhs_nrow[] */
137267602552SJunchao Zhang         mumps->nloc_rhs = 0;
13739566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRanges(A, &ranges));
137467602552SJunchao Zhang         for (j = 0; j < ompsize; j++) {
137567602552SJunchao Zhang           mumps->rhs_nrow[j] = ranges[petsc_ranks[j] + 1] - ranges[petsc_ranks[j]];
137667602552SJunchao Zhang           mumps->nloc_rhs += mumps->rhs_nrow[j];
137767602552SJunchao Zhang         }
13789566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(mumps->nloc_rhs, &mumps->irhs_loc));
137967602552SJunchao Zhang         for (j = k = 0; j < ompsize; j++) {
138067602552SJunchao 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 */
138167602552SJunchao Zhang         }
138267602552SJunchao Zhang 
13839566063dSJacob Faibussowitsch         PetscCall(PetscFree2(omp_ranks, petsc_ranks));
13849566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_free(&petsc_group));
13859566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_free(&omp_group));
138667602552SJunchao Zhang       }
138767602552SJunchao Zhang 
138867602552SJunchao Zhang       /* Realloc buffers when current nrhs is bigger than what we have met */
138967602552SJunchao Zhang       if (nrhs > mumps->max_nrhs) {
13909566063dSJacob Faibussowitsch         PetscCall(PetscFree2(mumps->rhs_loc, mumps->rhs_recvbuf));
13919566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(mumps->nloc_rhs * nrhs, &mumps->rhs_loc, mumps->nloc_rhs * nrhs, &mumps->rhs_recvbuf));
139267602552SJunchao Zhang         mumps->max_nrhs = nrhs;
139367602552SJunchao Zhang       }
139467602552SJunchao Zhang 
139567602552SJunchao Zhang       /* Setup recvcounts[], disps[], recvbuf on omp rank 0 for the upcoming MPI_Gatherv */
13969566063dSJacob Faibussowitsch       for (j = 0; j < ompsize; j++) PetscCall(PetscMPIIntCast(mumps->rhs_nrow[j] * nrhs, &mumps->rhs_recvcounts[j]));
139767602552SJunchao Zhang       mumps->rhs_disps[0] = 0;
139867602552SJunchao Zhang       for (j = 1; j < ompsize; j++) {
139967602552SJunchao Zhang         mumps->rhs_disps[j] = mumps->rhs_disps[j - 1] + mumps->rhs_recvcounts[j - 1];
140008401ef6SPierre Jolivet         PetscCheck(mumps->rhs_disps[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "PetscMPIInt overflow!");
140167602552SJunchao Zhang       }
140267602552SJunchao Zhang       recvbuf = (nrhs == 1) ? mumps->rhs_loc : mumps->rhs_recvbuf; /* Directly use rhs_loc[] as recvbuf. Single rhs is common in Ax=b */
140367602552SJunchao Zhang     }
140467602552SJunchao Zhang 
14059566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(m * nrhs, &sendcount));
14069566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(array, sendcount, MPIU_SCALAR, recvbuf, mumps->rhs_recvcounts, mumps->rhs_disps, MPIU_SCALAR, 0, mumps->omp_comm));
140767602552SJunchao Zhang 
140867602552SJunchao Zhang     if (mumps->is_omp_master) {
140967602552SJunchao Zhang       if (nrhs > 1) { /* Copy & re-arrange data from rhs_recvbuf[] to mumps->rhs_loc[] only when there are multiple rhs */
141067602552SJunchao Zhang         PetscScalar *dst, *dstbase = mumps->rhs_loc;
141167602552SJunchao Zhang         for (j = 0; j < ompsize; j++) {
141267602552SJunchao Zhang           const PetscScalar *src = mumps->rhs_recvbuf + mumps->rhs_disps[j];
141367602552SJunchao Zhang           dst                    = dstbase;
141467602552SJunchao Zhang           for (i = 0; i < nrhs; i++) {
14159566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(dst, src, mumps->rhs_nrow[j]));
141667602552SJunchao Zhang             src += mumps->rhs_nrow[j];
141767602552SJunchao Zhang             dst += mumps->nloc_rhs;
141867602552SJunchao Zhang           }
141967602552SJunchao Zhang           dstbase += mumps->rhs_nrow[j];
142067602552SJunchao Zhang         }
142167602552SJunchao Zhang       }
142267602552SJunchao Zhang       mumps->id.rhs_loc = (MumpsScalar *)mumps->rhs_loc;
142367602552SJunchao Zhang     }
142467602552SJunchao Zhang #endif /* PETSC_HAVE_OPENMP_SUPPORT */
142567602552SJunchao Zhang   }
14266497c311SBarry Smith   mumps->id.nrhs     = (PetscMUMPSInt)nrhs;
14276497c311SBarry Smith   mumps->id.nloc_rhs = (PetscMUMPSInt)mumps->nloc_rhs;
142867602552SJunchao Zhang   mumps->id.lrhs_loc = mumps->nloc_rhs;
142967602552SJunchao Zhang   mumps->id.irhs_loc = mumps->irhs_loc;
14303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143167602552SJunchao Zhang }
143267602552SJunchao Zhang 
143366976f2fSJacob Faibussowitsch static PetscErrorCode MatSolve_MUMPS(Mat A, Vec b, Vec x)
1434d71ae5a4SJacob Faibussowitsch {
1435e69c285eSBarry Smith   Mat_MUMPS         *mumps  = (Mat_MUMPS *)A->data;
143625aac85cSJunchao Zhang   const PetscScalar *rarray = NULL;
1437d54de34fSKris Buschelman   PetscScalar       *array;
1438329ec9b3SHong Zhang   IS                 is_iden, is_petsc;
1439329ec9b3SHong Zhang   PetscInt           i;
1440cc86f929SStefano Zampini   PetscBool          second_solve = PETSC_FALSE;
1441883f2eb9SBarry Smith   static PetscBool   cite1 = PETSC_FALSE, cite2 = PETSC_FALSE;
1442397b6df1SKris Buschelman 
1443397b6df1SKris Buschelman   PetscFunctionBegin;
14449371c9d4SSatish 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 "
14459371c9d4SSatish Balay                                    "Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",
14469371c9d4SSatish Balay                                    &cite1));
14479371c9d4SSatish 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 "
14489371c9d4SSatish Balay                                    "Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",
14499371c9d4SSatish Balay                                    &cite2));
14502aca8efcSHong Zhang 
1451f480ea8aSBarry Smith   PetscCall(VecFlag(x, A->factorerrortype));
1452603e8f96SBarry Smith   if (A->factorerrortype) {
14539566063dSJacob 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)));
14543ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
14552aca8efcSHong Zhang   }
14562aca8efcSHong Zhang 
1457a5e57a09SHong Zhang   mumps->id.nrhs = 1;
14582d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
145925aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
146067602552SJunchao Zhang       mumps->id.ICNTL(20) = 10; /* dense distributed RHS */
14619566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(b, &rarray));
14629566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetUpDistRHSInfo(A, 1, rarray));
146325aac85cSJunchao Zhang     } else {
146441ffd417SStefano Zampini       mumps->id.ICNTL(20) = 0; /* dense centralized RHS; Scatter b into a sequential rhs vector*/
14659566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mumps->scat_rhs, b, mumps->b_seq, INSERT_VALUES, SCATTER_FORWARD));
14669566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mumps->scat_rhs, b, mumps->b_seq, INSERT_VALUES, SCATTER_FORWARD));
146767602552SJunchao Zhang       if (!mumps->myid) {
14689566063dSJacob Faibussowitsch         PetscCall(VecGetArray(mumps->b_seq, &array));
146967602552SJunchao Zhang         mumps->id.rhs = (MumpsScalar *)array;
147067602552SJunchao Zhang       }
147125aac85cSJunchao Zhang     }
14723ab56b82SJunchao Zhang   } else {                   /* petsc_size == 1 */
147367602552SJunchao Zhang     mumps->id.ICNTL(20) = 0; /* dense centralized RHS */
14749566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, x));
14759566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &array));
1476940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar *)array;
1477397b6df1SKris Buschelman   }
1478397b6df1SKris Buschelman 
1479cc86f929SStefano Zampini   /*
1480cc86f929SStefano Zampini      handle condensation step of Schur complement (if any)
1481cc86f929SStefano Zampini      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
1482cc86f929SStefano Zampini      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
1483cc86f929SStefano Zampini      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
1484cc86f929SStefano Zampini      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
1485cc86f929SStefano Zampini   */
14863e5b40d0SPierre Jolivet   if (mumps->id.size_schur > 0) {
148708401ef6SPierre Jolivet     PetscCheck(mumps->petsc_size <= 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Parallel Schur complements not yet supported from PETSc");
14883e5b40d0SPierre Jolivet     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
1489cc86f929SStefano Zampini       second_solve = PETSC_TRUE;
14909566063dSJacob Faibussowitsch       PetscCall(MatMumpsHandleSchur_Private(A, PETSC_FALSE));
14913e5b40d0SPierre Jolivet       mumps->id.ICNTL(26) = 1; /* condensation phase */
14923e5b40d0SPierre Jolivet     } else if (mumps->id.ICNTL(26) == 1) PetscCall(MatMumpsHandleSchur_Private(A, PETSC_FALSE));
1493cc86f929SStefano Zampini   }
1494397b6df1SKris Buschelman   /* solve phase */
1495a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
14963ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
14979261f6e4SBarry Smith   PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in solve: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
1498397b6df1SKris Buschelman 
1499b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
15001baa6e33SBarry Smith   if (second_solve) PetscCall(MatMumpsHandleSchur_Private(A, PETSC_TRUE));
15013e5b40d0SPierre Jolivet   else if (mumps->id.ICNTL(26) == 1) {
15023e5b40d0SPierre Jolivet     PetscCall(MatMumpsSolveSchur_Private(A));
15033e5b40d0SPierre Jolivet     for (i = 0; i < mumps->id.size_schur; ++i) {
15043e5b40d0SPierre Jolivet #if !defined(PETSC_USE_COMPLEX)
15053e5b40d0SPierre Jolivet       PetscScalar val = mumps->id.redrhs[i];
15063e5b40d0SPierre Jolivet #else
15073e5b40d0SPierre Jolivet       PetscScalar val = mumps->id.redrhs[i].r + PETSC_i * mumps->id.redrhs[i].i;
15083e5b40d0SPierre Jolivet #endif
15093e5b40d0SPierre Jolivet       array[mumps->id.listvar_schur[i] - 1] = val;
15103e5b40d0SPierre Jolivet     }
15113e5b40d0SPierre Jolivet   }
1512b5fa320bSStefano Zampini 
15132d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) { /* convert mumps distributed solution to petsc mpi x */
1514a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
1515a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
15169566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&mumps->scat_sol));
1517397b6df1SKris Buschelman     }
1518a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
1519a6053eceSJunchao Zhang       PetscInt *isol2_loc = NULL;
15209566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, mumps->id.lsol_loc, 0, 1, &is_iden)); /* from */
15219566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mumps->id.lsol_loc, &isol2_loc));
1522a6053eceSJunchao 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 */
15239566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, mumps->id.lsol_loc, isol2_loc, PETSC_OWN_POINTER, &is_petsc)); /* to */
15249566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(mumps->x_seq, is_iden, x, is_petsc, &mumps->scat_sol));
15259566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_iden));
15269566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_petsc));
1527a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
1528397b6df1SKris Buschelman     }
1529a5e57a09SHong Zhang 
15309566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mumps->scat_sol, mumps->x_seq, x, INSERT_VALUES, SCATTER_FORWARD));
15319566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mumps->scat_sol, mumps->x_seq, x, INSERT_VALUES, SCATTER_FORWARD));
1532329ec9b3SHong Zhang   }
1533353d7d71SJunchao Zhang 
153467602552SJunchao Zhang   if (mumps->petsc_size > 1) {
153525aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
15369566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(b, &rarray));
153725aac85cSJunchao Zhang     } else if (!mumps->myid) {
15389566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mumps->b_seq, &array));
153925aac85cSJunchao Zhang     }
15409566063dSJacob Faibussowitsch   } else PetscCall(VecRestoreArray(x, &array));
1541353d7d71SJunchao Zhang 
154264412097SPierre Jolivet   PetscCall(PetscLogFlops(2.0 * PetscMax(0, (mumps->id.INFO(28) >= 0 ? mumps->id.INFO(28) : -1000000 * mumps->id.INFO(28)) - A->cmap->n)));
15433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1544397b6df1SKris Buschelman }
1545397b6df1SKris Buschelman 
154666976f2fSJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_MUMPS(Mat A, Vec b, Vec x)
1547d71ae5a4SJacob Faibussowitsch {
1548e69c285eSBarry Smith   Mat_MUMPS          *mumps = (Mat_MUMPS *)A->data;
1549338d3105SPierre Jolivet   const PetscMUMPSInt value = mumps->id.ICNTL(9);
155051d5961aSHong Zhang 
155151d5961aSHong Zhang   PetscFunctionBegin;
1552a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
15539566063dSJacob Faibussowitsch   PetscCall(MatSolve_MUMPS(A, b, x));
1554338d3105SPierre Jolivet   mumps->id.ICNTL(9) = value;
15553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155651d5961aSHong Zhang }
155751d5961aSHong Zhang 
155866976f2fSJacob Faibussowitsch static PetscErrorCode MatMatSolve_MUMPS(Mat A, Mat B, Mat X)
1559d71ae5a4SJacob Faibussowitsch {
1560b8491c3eSStefano Zampini   Mat                Bt = NULL;
1561a6053eceSJunchao Zhang   PetscBool          denseX, denseB, flg, flgT;
1562e69c285eSBarry Smith   Mat_MUMPS         *mumps = (Mat_MUMPS *)A->data;
1563917c3dccSPierre Jolivet   PetscInt           i, nrhs, M, nrhsM;
15641683a169SBarry Smith   PetscScalar       *array;
15651683a169SBarry Smith   const PetscScalar *rbray;
1566a6053eceSJunchao Zhang   PetscInt           lsol_loc, nlsol_loc, *idxx, iidx = 0;
1567a6053eceSJunchao Zhang   PetscMUMPSInt     *isol_loc, *isol_loc_save;
15681683a169SBarry Smith   PetscScalar       *bray, *sol_loc, *sol_loc_save;
1569be818407SHong Zhang   IS                 is_to, is_from;
1570beae5ec0SHong Zhang   PetscInt           k, proc, j, m, myrstart;
1571be818407SHong Zhang   const PetscInt    *rstart;
157267602552SJunchao Zhang   Vec                v_mpi, msol_loc;
157367602552SJunchao Zhang   VecScatter         scat_sol;
157467602552SJunchao Zhang   Vec                b_seq;
157567602552SJunchao Zhang   VecScatter         scat_rhs;
1576be818407SHong Zhang   PetscScalar       *aa;
1577be818407SHong Zhang   PetscInt           spnr, *ia, *ja;
1578d56c302dSHong Zhang   Mat_MPIAIJ        *b = NULL;
1579bda8bf91SBarry Smith 
1580e0b74bf9SHong Zhang   PetscFunctionBegin;
15819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &denseX, MATSEQDENSE, MATMPIDENSE, NULL));
158228b400f6SJacob Faibussowitsch   PetscCheck(denseX, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix");
1583be818407SHong Zhang 
15849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &denseB, MATSEQDENSE, MATMPIDENSE, NULL));
1585a6053eceSJunchao Zhang   if (denseB) {
158608401ef6SPierre Jolivet     PetscCheck(B->rmap->n == X->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix B and X must have same row distribution");
1587be818407SHong Zhang     mumps->id.ICNTL(20) = 0; /* dense RHS */
15880e6b8875SHong Zhang   } else {                   /* sparse B */
158908401ef6SPierre Jolivet     PetscCheck(X != B, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_IDN, "X and B must be different matrices");
1590013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &flgT));
1591da81f932SPierre Jolivet     if (flgT) { /* input B is transpose of actual RHS matrix,
15920e6b8875SHong Zhang                  because mumps requires sparse compressed COLUMN storage! See MatMatTransposeSolve_MUMPS() */
15939566063dSJacob Faibussowitsch       PetscCall(MatTransposeGetMat(B, &Bt));
1594013e2dc7SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONG, "Matrix B must be MATTRANSPOSEVIRTUAL matrix");
1595be818407SHong Zhang     mumps->id.ICNTL(20) = 1; /* sparse RHS */
1596b8491c3eSStefano Zampini   }
159787b22cf4SHong Zhang 
15989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, &nrhs));
1599917c3dccSPierre Jolivet   PetscCall(PetscIntMultError(nrhs, M, &nrhsM));
16006497c311SBarry Smith   mumps->id.nrhs = (PetscMUMPSInt)nrhs;
16016497c311SBarry Smith   mumps->id.lrhs = (PetscMUMPSInt)M;
16022b691707SHong Zhang   mumps->id.rhs  = NULL;
16039481e6e9SHong Zhang 
16042d4298aeSJunchao Zhang   if (mumps->petsc_size == 1) {
1605b8491c3eSStefano Zampini     PetscScalar *aa;
1606b8491c3eSStefano Zampini     PetscInt     spnr, *ia, *ja;
1607e94cce23SStefano Zampini     PetscBool    second_solve = PETSC_FALSE;
1608b8491c3eSStefano Zampini 
16099566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(X, &array));
1610b8491c3eSStefano Zampini     mumps->id.rhs = (MumpsScalar *)array;
16112b691707SHong Zhang 
1612a6053eceSJunchao Zhang     if (denseB) {
16132b691707SHong Zhang       /* copy B to X */
16149566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(B, &rbray));
1615917c3dccSPierre Jolivet       PetscCall(PetscArraycpy(array, rbray, nrhsM));
16169566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(B, &rbray));
16172b691707SHong Zhang     } else { /* sparse B */
16189566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(Bt, &aa));
16199566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Bt, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
162028b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
16219566063dSJacob Faibussowitsch       PetscCall(PetscMUMPSIntCSRCast(mumps, spnr, ia, ja, &mumps->id.irhs_ptr, &mumps->id.irhs_sparse, &mumps->id.nz_rhs));
1622b8491c3eSStefano Zampini       mumps->id.rhs_sparse = (MumpsScalar *)aa;
1623b8491c3eSStefano Zampini     }
1624e94cce23SStefano Zampini     /* handle condensation step of Schur complement (if any) */
16253e5b40d0SPierre Jolivet     if (mumps->id.size_schur > 0) {
16263e5b40d0SPierre Jolivet       if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
1627e94cce23SStefano Zampini         second_solve = PETSC_TRUE;
16289566063dSJacob Faibussowitsch         PetscCall(MatMumpsHandleSchur_Private(A, PETSC_FALSE));
16293e5b40d0SPierre Jolivet         mumps->id.ICNTL(26) = 1; /* condensation phase */
16303e5b40d0SPierre Jolivet       } else if (mumps->id.ICNTL(26) == 1) PetscCall(MatMumpsHandleSchur_Private(A, PETSC_FALSE));
1631e94cce23SStefano Zampini     }
16322cd7d884SHong Zhang     /* solve phase */
16332cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
16343ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
16359261f6e4SBarry Smith     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in solve: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
1636b5fa320bSStefano Zampini 
1637b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
16381baa6e33SBarry Smith     if (second_solve) PetscCall(MatMumpsHandleSchur_Private(A, PETSC_TRUE));
16393e5b40d0SPierre Jolivet     else if (mumps->id.ICNTL(26) == 1) {
16403e5b40d0SPierre Jolivet       PetscCall(MatMumpsSolveSchur_Private(A));
16413e5b40d0SPierre Jolivet       for (j = 0; j < nrhs; ++j)
16423e5b40d0SPierre Jolivet         for (i = 0; i < mumps->id.size_schur; ++i) {
16433e5b40d0SPierre Jolivet #if !defined(PETSC_USE_COMPLEX)
16443e5b40d0SPierre Jolivet           PetscScalar val = mumps->id.redrhs[i + j * mumps->id.lredrhs];
16453e5b40d0SPierre Jolivet #else
16463e5b40d0SPierre Jolivet           PetscScalar val = mumps->id.redrhs[i + j * mumps->id.lredrhs].r + PETSC_i * mumps->id.redrhs[i + j * mumps->id.lredrhs].i;
16473e5b40d0SPierre Jolivet #endif
16483e5b40d0SPierre Jolivet           array[mumps->id.listvar_schur[i] - 1 + j * M] = val;
16493e5b40d0SPierre Jolivet         }
16503e5b40d0SPierre Jolivet     }
1651a6053eceSJunchao Zhang     if (!denseB) { /* sparse B */
16529566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(Bt, &aa));
16539566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Bt, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
165428b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot restore IJ structure");
1655b8491c3eSStefano Zampini     }
16569566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X, &array));
16573ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1658be818407SHong Zhang   }
1659801fbe65SHong Zhang 
16602ef1f0ffSBarry Smith   /* parallel case: MUMPS requires rhs B to be centralized on the host! */
166150a7cd33SPierre Jolivet   PetscCheck(!mumps->id.ICNTL(19), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Parallel Schur complements not yet supported from PETSc");
1662241dbb5eSStefano Zampini 
1663beae5ec0SHong Zhang   /* create msol_loc to hold mumps local solution */
16641683a169SBarry Smith   isol_loc_save = mumps->id.isol_loc; /* save it for MatSolve() */
16651683a169SBarry Smith   sol_loc_save  = (PetscScalar *)mumps->id.sol_loc;
1666801fbe65SHong Zhang 
1667a1dfcbd9SJunchao Zhang   lsol_loc = mumps->id.lsol_loc;
1668917c3dccSPierre Jolivet   PetscCall(PetscIntMultError(nrhs, lsol_loc, &nlsol_loc)); /* length of sol_loc */
16699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nlsol_loc, &sol_loc, lsol_loc, &isol_loc));
1670940cd9d6SSatish Balay   mumps->id.sol_loc  = (MumpsScalar *)sol_loc;
1671801fbe65SHong Zhang   mumps->id.isol_loc = isol_loc;
1672801fbe65SHong Zhang 
16739566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nlsol_loc, (PetscScalar *)sol_loc, &msol_loc));
16742cd7d884SHong Zhang 
167567602552SJunchao Zhang   if (denseB) {
167625aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
167767602552SJunchao Zhang       mumps->id.ICNTL(20) = 10; /* dense distributed RHS */
16789566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(B, &rbray));
16799566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetUpDistRHSInfo(A, nrhs, rbray));
16809566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(B, &rbray));
16819566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, &m, NULL));
1682917c3dccSPierre Jolivet       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)B), 1, nrhs * m, nrhsM, NULL, &v_mpi));
168325aac85cSJunchao Zhang     } else {
168425aac85cSJunchao Zhang       mumps->id.ICNTL(20) = 0; /* dense centralized RHS */
168580577c12SJunchao Zhang       /* TODO: Because of non-contiguous indices, the created vecscatter scat_rhs is not done in MPI_Gather, resulting in
168680577c12SJunchao Zhang         very inefficient communication. An optimization is to use VecScatterCreateToZero to gather B to rank 0. Then on rank
168780577c12SJunchao Zhang         0, re-arrange B into desired order, which is a local operation.
168880577c12SJunchao Zhang       */
168980577c12SJunchao Zhang 
169067602552SJunchao Zhang       /* scatter v_mpi to b_seq because MUMPS before 5.3.0 only supports centralized rhs */
1691be818407SHong Zhang       /* wrap dense rhs matrix B into a vector v_mpi */
16929566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, &m, NULL));
16939566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArray(B, &bray));
1694917c3dccSPierre Jolivet       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)B), 1, nrhs * m, nrhsM, (const PetscScalar *)bray, &v_mpi));
16959566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArray(B, &bray));
16962b691707SHong Zhang 
1697be818407SHong Zhang       /* scatter v_mpi to b_seq in proc[0]. MUMPS requires rhs to be centralized on the host! */
1698801fbe65SHong Zhang       if (!mumps->myid) {
1699beae5ec0SHong Zhang         PetscInt *idx;
1700beae5ec0SHong Zhang         /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B */
1701917c3dccSPierre Jolivet         PetscCall(PetscMalloc1(nrhsM, &idx));
17029566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRanges(B, &rstart));
1703917c3dccSPierre Jolivet         for (proc = 0, k = 0; proc < mumps->petsc_size; proc++) {
1704be818407SHong Zhang           for (j = 0; j < nrhs; j++) {
1705beae5ec0SHong Zhang             for (i = rstart[proc]; i < rstart[proc + 1]; i++) idx[k++] = j * M + i;
1706be818407SHong Zhang           }
1707be818407SHong Zhang         }
1708be818407SHong Zhang 
1709917c3dccSPierre Jolivet         PetscCall(VecCreateSeq(PETSC_COMM_SELF, nrhsM, &b_seq));
1710917c3dccSPierre Jolivet         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrhsM, idx, PETSC_OWN_POINTER, &is_to));
1711917c3dccSPierre Jolivet         PetscCall(ISCreateStride(PETSC_COMM_SELF, nrhsM, 0, 1, &is_from));
1712801fbe65SHong Zhang       } else {
17139566063dSJacob Faibussowitsch         PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &b_seq));
17149566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is_to));
17159566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is_from));
1716801fbe65SHong Zhang       }
17179566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(v_mpi, is_from, b_seq, is_to, &scat_rhs));
17189566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(scat_rhs, v_mpi, b_seq, INSERT_VALUES, SCATTER_FORWARD));
17199566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_to));
17209566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_from));
17219566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(scat_rhs, v_mpi, b_seq, INSERT_VALUES, SCATTER_FORWARD));
1722801fbe65SHong Zhang 
1723801fbe65SHong Zhang       if (!mumps->myid) { /* define rhs on the host */
17249566063dSJacob Faibussowitsch         PetscCall(VecGetArray(b_seq, &bray));
1725940cd9d6SSatish Balay         mumps->id.rhs = (MumpsScalar *)bray;
17269566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(b_seq, &bray));
1727801fbe65SHong Zhang       }
172825aac85cSJunchao Zhang     }
17292b691707SHong Zhang   } else { /* sparse B */
17302b691707SHong Zhang     b = (Mat_MPIAIJ *)Bt->data;
17312b691707SHong Zhang 
1732be818407SHong Zhang     /* wrap dense X into a vector v_mpi */
17339566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(X, &m, NULL));
17349566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(X, &bray));
1735917c3dccSPierre Jolivet     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)X), 1, nrhs * m, nrhsM, (const PetscScalar *)bray, &v_mpi));
17369566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X, &bray));
17372b691707SHong Zhang 
17382b691707SHong Zhang     if (!mumps->myid) {
17399566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(b->A, &aa));
17409566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(b->A, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
174128b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
17429566063dSJacob Faibussowitsch       PetscCall(PetscMUMPSIntCSRCast(mumps, spnr, ia, ja, &mumps->id.irhs_ptr, &mumps->id.irhs_sparse, &mumps->id.nz_rhs));
17432b691707SHong Zhang       mumps->id.rhs_sparse = (MumpsScalar *)aa;
17442b691707SHong Zhang     } else {
17452b691707SHong Zhang       mumps->id.irhs_ptr    = NULL;
17462b691707SHong Zhang       mumps->id.irhs_sparse = NULL;
17472b691707SHong Zhang       mumps->id.nz_rhs      = 0;
17482b691707SHong Zhang       mumps->id.rhs_sparse  = NULL;
17492b691707SHong Zhang     }
17502b691707SHong Zhang   }
17512b691707SHong Zhang 
1752801fbe65SHong Zhang   /* solve phase */
1753801fbe65SHong Zhang   mumps->id.job = JOB_SOLVE;
17543ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
17559261f6e4SBarry Smith   PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in solve: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
1756801fbe65SHong Zhang 
1757334c5f61SHong Zhang   /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
17589566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(X, &array));
17599566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(v_mpi, array));
1760801fbe65SHong Zhang 
1761334c5f61SHong Zhang   /* create scatter scat_sol */
17629566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(X, &rstart));
1763beae5ec0SHong Zhang   /* iidx: index for scatter mumps solution to petsc X */
1764beae5ec0SHong Zhang 
17659566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, nlsol_loc, 0, 1, &is_from));
17669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nlsol_loc, &idxx));
1767beae5ec0SHong Zhang   for (i = 0; i < lsol_loc; i++) {
1768beae5ec0SHong 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 */
1769beae5ec0SHong Zhang 
17702d4298aeSJunchao Zhang     for (proc = 0; proc < mumps->petsc_size; proc++) {
1771beae5ec0SHong Zhang       if (isol_loc[i] >= rstart[proc] && isol_loc[i] < rstart[proc + 1]) {
1772beae5ec0SHong Zhang         myrstart = rstart[proc];
1773beae5ec0SHong Zhang         k        = isol_loc[i] - myrstart;          /* local index on 1st column of petsc vector X */
1774beae5ec0SHong Zhang         iidx     = k + myrstart * nrhs;             /* maps mumps isol_loc[i] to petsc index in X */
1775beae5ec0SHong Zhang         m        = rstart[proc + 1] - rstart[proc]; /* rows of X for this proc */
1776beae5ec0SHong Zhang         break;
1777be818407SHong Zhang       }
1778be818407SHong Zhang     }
1779be818407SHong Zhang 
1780beae5ec0SHong Zhang     for (j = 0; j < nrhs; j++) idxx[i + j * lsol_loc] = iidx + j * m;
1781801fbe65SHong Zhang   }
17829566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nlsol_loc, idxx, PETSC_COPY_VALUES, &is_to));
17839566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(msol_loc, is_from, v_mpi, is_to, &scat_sol));
17849566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat_sol, msol_loc, v_mpi, INSERT_VALUES, SCATTER_FORWARD));
17859566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_from));
17869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_to));
17879566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat_sol, msol_loc, v_mpi, INSERT_VALUES, SCATTER_FORWARD));
17889566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(X, &array));
178971aed81dSHong Zhang 
179071aed81dSHong Zhang   /* free spaces */
17911683a169SBarry Smith   mumps->id.sol_loc  = (MumpsScalar *)sol_loc_save;
179271aed81dSHong Zhang   mumps->id.isol_loc = isol_loc_save;
179371aed81dSHong Zhang 
17949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sol_loc, isol_loc));
17959566063dSJacob Faibussowitsch   PetscCall(PetscFree(idxx));
17969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&msol_loc));
17979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v_mpi));
1798a6053eceSJunchao Zhang   if (!denseB) {
17992b691707SHong Zhang     if (!mumps->myid) {
1800d56c302dSHong Zhang       b = (Mat_MPIAIJ *)Bt->data;
18019566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(b->A, &aa));
18029566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(b->A, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
180328b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot restore IJ structure");
18042b691707SHong Zhang     }
18052b691707SHong Zhang   } else {
180625aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) {
18079566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b_seq));
18089566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&scat_rhs));
180925aac85cSJunchao Zhang     }
18102b691707SHong Zhang   }
18119566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scat_sol));
181257508eceSPierre Jolivet   PetscCall(PetscLogFlops(nrhs * PetscMax(0, 2.0 * (mumps->id.INFO(28) >= 0 ? mumps->id.INFO(28) : -1000000 * mumps->id.INFO(28)) - A->cmap->n)));
18133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1814e0b74bf9SHong Zhang }
1815e0b74bf9SHong Zhang 
181666976f2fSJacob Faibussowitsch static PetscErrorCode MatMatSolveTranspose_MUMPS(Mat A, Mat B, Mat X)
1817d71ae5a4SJacob Faibussowitsch {
1818b18964edSHong Zhang   Mat_MUMPS          *mumps = (Mat_MUMPS *)A->data;
1819338d3105SPierre Jolivet   const PetscMUMPSInt value = mumps->id.ICNTL(9);
1820b18964edSHong Zhang 
1821b18964edSHong Zhang   PetscFunctionBegin;
1822b18964edSHong Zhang   mumps->id.ICNTL(9) = 0;
1823b18964edSHong Zhang   PetscCall(MatMatSolve_MUMPS(A, B, X));
1824338d3105SPierre Jolivet   mumps->id.ICNTL(9) = value;
18253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1826b18964edSHong Zhang }
1827b18964edSHong Zhang 
182866976f2fSJacob Faibussowitsch static PetscErrorCode MatMatTransposeSolve_MUMPS(Mat A, Mat Bt, Mat X)
1829d71ae5a4SJacob Faibussowitsch {
1830eb3ef3b2SHong Zhang   PetscBool flg;
1831eb3ef3b2SHong Zhang   Mat       B;
1832eb3ef3b2SHong Zhang 
1833eb3ef3b2SHong Zhang   PetscFunctionBegin;
18349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)Bt, &flg, MATSEQAIJ, MATMPIAIJ, NULL));
183528b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)Bt), PETSC_ERR_ARG_WRONG, "Matrix Bt must be MATAIJ matrix");
1836eb3ef3b2SHong Zhang 
1837eb3ef3b2SHong Zhang   /* Create B=Bt^T that uses Bt's data structure */
18389566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(Bt, &B));
1839eb3ef3b2SHong Zhang 
18409566063dSJacob Faibussowitsch   PetscCall(MatMatSolve_MUMPS(A, B, X));
18419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
18423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1843eb3ef3b2SHong Zhang }
1844eb3ef3b2SHong Zhang 
1845ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1846a58c3f20SHong Zhang /*
1847a58c3f20SHong Zhang   input:
1848a58c3f20SHong Zhang    F:        numeric factor
1849a58c3f20SHong Zhang   output:
1850a58c3f20SHong Zhang    nneg:     total number of negative pivots
185119d49a3bSHong Zhang    nzero:    total number of zero pivots
185219d49a3bSHong Zhang    npos:     (global dimension of F) - nneg - nzero
1853a58c3f20SHong Zhang */
185466976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F, PetscInt *nneg, PetscInt *nzero, PetscInt *npos)
1855d71ae5a4SJacob Faibussowitsch {
1856e69c285eSBarry Smith   Mat_MUMPS  *mumps = (Mat_MUMPS *)F->data;
1857c1490034SHong Zhang   PetscMPIInt size;
1858a58c3f20SHong Zhang 
1859a58c3f20SHong Zhang   PetscFunctionBegin;
18609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)F), &size));
1861bcb30aebSHong 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 */
1862aed4548fSBarry 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));
1863ed85ac9fSHong Zhang 
1864710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1865ed85ac9fSHong Zhang   if (nzero || npos) {
186608401ef6SPierre 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");
1867710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1868710ac8efSHong Zhang     if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1869a58c3f20SHong Zhang   }
18703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1871a58c3f20SHong Zhang }
187219d49a3bSHong Zhang #endif
1873a58c3f20SHong Zhang 
187466976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGatherNonzerosOnMaster(MatReuse reuse, Mat_MUMPS *mumps)
1875d71ae5a4SJacob Faibussowitsch {
18766497c311SBarry Smith   PetscMPIInt    nreqs;
1877a6053eceSJunchao Zhang   PetscMUMPSInt *irn, *jcn;
1878a6053eceSJunchao Zhang   PetscMPIInt    count;
18796497c311SBarry Smith   PetscCount     totnnz, remain;
1880a6053eceSJunchao Zhang   const PetscInt osize = mumps->omp_comm_size;
1881a6053eceSJunchao Zhang   PetscScalar   *val;
18823ab56b82SJunchao Zhang 
18833ab56b82SJunchao Zhang   PetscFunctionBegin;
1884a6053eceSJunchao Zhang   if (osize > 1) {
18853ab56b82SJunchao Zhang     if (reuse == MAT_INITIAL_MATRIX) {
18863ab56b82SJunchao Zhang       /* master first gathers counts of nonzeros to receive */
18879566063dSJacob Faibussowitsch       if (mumps->is_omp_master) PetscCall(PetscMalloc1(osize, &mumps->recvcount));
18889566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gather(&mumps->nnz, 1, MPIU_INT64, mumps->recvcount, 1, MPIU_INT64, 0 /*master*/, mumps->omp_comm));
18893ab56b82SJunchao Zhang 
1890a6053eceSJunchao Zhang       /* Then each computes number of send/recvs */
18913ab56b82SJunchao Zhang       if (mumps->is_omp_master) {
1892a6053eceSJunchao Zhang         /* Start from 1 since self communication is not done in MPI */
1893a6053eceSJunchao Zhang         nreqs = 0;
18946497c311SBarry Smith         for (PetscMPIInt i = 1; i < osize; i++) nreqs += (mumps->recvcount[i] + PETSC_MPI_INT_MAX - 1) / PETSC_MPI_INT_MAX;
1895a6053eceSJunchao Zhang       } else {
18966497c311SBarry Smith         nreqs = (PetscMPIInt)(((mumps->nnz + PETSC_MPI_INT_MAX - 1) / PETSC_MPI_INT_MAX));
18973ab56b82SJunchao Zhang       }
189835cb6cd3SPierre Jolivet       PetscCall(PetscMalloc1(nreqs * 3, &mumps->reqs)); /* Triple the requests since we send irn, jcn and val separately */
18993ab56b82SJunchao Zhang 
1900a6053eceSJunchao Zhang       /* The following code is doing a very simple thing: omp_master rank gathers irn/jcn/val from others.
1901a6053eceSJunchao Zhang          MPI_Gatherv would be enough if it supports big counts > 2^31-1. Since it does not, and mumps->nnz
1902a6053eceSJunchao Zhang          might be a prime number > 2^31-1, we have to slice the message. Note omp_comm_size
1903a6053eceSJunchao Zhang          is very small, the current approach should have no extra overhead compared to MPI_Gatherv.
1904a6053eceSJunchao Zhang        */
1905a6053eceSJunchao Zhang       nreqs = 0; /* counter for actual send/recvs */
19063ab56b82SJunchao Zhang       if (mumps->is_omp_master) {
19076497c311SBarry Smith         totnnz = 0;
19086497c311SBarry Smith 
19096497c311SBarry Smith         for (PetscMPIInt i = 0; i < osize; i++) totnnz += mumps->recvcount[i]; /* totnnz = sum of nnz over omp_comm */
19109566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(totnnz, &irn, totnnz, &jcn));
19119566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(totnnz, &val));
1912a6053eceSJunchao Zhang 
1913a6053eceSJunchao Zhang         /* Self communication */
19149566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(irn, mumps->irn, mumps->nnz));
19159566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jcn, mumps->jcn, mumps->nnz));
19169566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(val, mumps->val, mumps->nnz));
1917a6053eceSJunchao Zhang 
1918a6053eceSJunchao Zhang         /* Replace mumps->irn/jcn etc on master with the newly allocated bigger arrays */
19199566063dSJacob Faibussowitsch         PetscCall(PetscFree2(mumps->irn, mumps->jcn));
19209566063dSJacob Faibussowitsch         PetscCall(PetscFree(mumps->val_alloc));
1921a6053eceSJunchao Zhang         mumps->nnz = totnnz;
19223ab56b82SJunchao Zhang         mumps->irn = irn;
19233ab56b82SJunchao Zhang         mumps->jcn = jcn;
1924a6053eceSJunchao Zhang         mumps->val = mumps->val_alloc = val;
1925a6053eceSJunchao Zhang 
1926a6053eceSJunchao Zhang         irn += mumps->recvcount[0]; /* recvcount[0] is old mumps->nnz on omp rank 0 */
1927a6053eceSJunchao Zhang         jcn += mumps->recvcount[0];
1928a6053eceSJunchao Zhang         val += mumps->recvcount[0];
1929a6053eceSJunchao Zhang 
1930a6053eceSJunchao Zhang         /* Remote communication */
19316497c311SBarry Smith         for (PetscMPIInt i = 1; i < osize; i++) {
19326497c311SBarry Smith           count  = (PetscMPIInt)PetscMin(mumps->recvcount[i], (PetscMPIInt)PETSC_MPI_INT_MAX);
1933a6053eceSJunchao Zhang           remain = mumps->recvcount[i] - count;
1934a6053eceSJunchao Zhang           while (count > 0) {
19356497c311SBarry Smith             PetscCallMPI(MPIU_Irecv(irn, count, MPIU_MUMPSINT, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
19366497c311SBarry Smith             PetscCallMPI(MPIU_Irecv(jcn, count, MPIU_MUMPSINT, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
19376497c311SBarry Smith             PetscCallMPI(MPIU_Irecv(val, count, MPIU_SCALAR, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1938a6053eceSJunchao Zhang             irn += count;
1939a6053eceSJunchao Zhang             jcn += count;
1940a6053eceSJunchao Zhang             val += count;
19416497c311SBarry Smith             count = (PetscMPIInt)PetscMin(remain, (PetscMPIInt)PETSC_MPI_INT_MAX);
1942a6053eceSJunchao Zhang             remain -= count;
1943a6053eceSJunchao Zhang           }
19443ab56b82SJunchao Zhang         }
19453ab56b82SJunchao Zhang       } else {
1946a6053eceSJunchao Zhang         irn    = mumps->irn;
1947a6053eceSJunchao Zhang         jcn    = mumps->jcn;
1948a6053eceSJunchao Zhang         val    = mumps->val;
19496497c311SBarry Smith         count  = (PetscMPIInt)PetscMin(mumps->nnz, (PetscMPIInt)PETSC_MPI_INT_MAX);
1950a6053eceSJunchao Zhang         remain = mumps->nnz - count;
1951a6053eceSJunchao Zhang         while (count > 0) {
19526497c311SBarry Smith           PetscCallMPI(MPIU_Isend(irn, count, MPIU_MUMPSINT, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
19536497c311SBarry Smith           PetscCallMPI(MPIU_Isend(jcn, count, MPIU_MUMPSINT, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
19546497c311SBarry Smith           PetscCallMPI(MPIU_Isend(val, count, MPIU_SCALAR, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1955a6053eceSJunchao Zhang           irn += count;
1956a6053eceSJunchao Zhang           jcn += count;
1957a6053eceSJunchao Zhang           val += count;
19586497c311SBarry Smith           count = (PetscMPIInt)PetscMin(remain, (PetscMPIInt)PETSC_MPI_INT_MAX);
1959a6053eceSJunchao Zhang           remain -= count;
19603ab56b82SJunchao Zhang         }
19613ab56b82SJunchao Zhang       }
1962a6053eceSJunchao Zhang     } else {
1963a6053eceSJunchao Zhang       nreqs = 0;
1964a6053eceSJunchao Zhang       if (mumps->is_omp_master) {
1965a6053eceSJunchao Zhang         val = mumps->val + mumps->recvcount[0];
19666497c311SBarry Smith         for (PetscMPIInt i = 1; i < osize; i++) { /* Remote communication only since self data is already in place */
19676497c311SBarry Smith           count  = (PetscMPIInt)PetscMin(mumps->recvcount[i], (PetscMPIInt)PETSC_MPI_INT_MAX);
1968a6053eceSJunchao Zhang           remain = mumps->recvcount[i] - count;
1969a6053eceSJunchao Zhang           while (count > 0) {
19706497c311SBarry Smith             PetscCallMPI(MPIU_Irecv(val, count, MPIU_SCALAR, i, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1971a6053eceSJunchao Zhang             val += count;
19726497c311SBarry Smith             count = (PetscMPIInt)PetscMin(remain, (PetscMPIInt)PETSC_MPI_INT_MAX);
1973a6053eceSJunchao Zhang             remain -= count;
1974a6053eceSJunchao Zhang           }
1975a6053eceSJunchao Zhang         }
1976a6053eceSJunchao Zhang       } else {
1977a6053eceSJunchao Zhang         val    = mumps->val;
19786497c311SBarry Smith         count  = (PetscMPIInt)PetscMin(mumps->nnz, (PetscMPIInt)PETSC_MPI_INT_MAX);
1979a6053eceSJunchao Zhang         remain = mumps->nnz - count;
1980a6053eceSJunchao Zhang         while (count > 0) {
19816497c311SBarry Smith           PetscCallMPI(MPIU_Isend(val, count, MPIU_SCALAR, 0, mumps->tag, mumps->omp_comm, &mumps->reqs[nreqs++]));
1982a6053eceSJunchao Zhang           val += count;
19836497c311SBarry Smith           count = (PetscMPIInt)PetscMin(remain, (PetscMPIInt)PETSC_MPI_INT_MAX);
1984a6053eceSJunchao Zhang           remain -= count;
1985a6053eceSJunchao Zhang         }
1986a6053eceSJunchao Zhang       }
1987a6053eceSJunchao Zhang     }
19889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nreqs, mumps->reqs, MPI_STATUSES_IGNORE));
1989a6053eceSJunchao Zhang     mumps->tag++; /* It is totally fine for above send/recvs to share one mpi tag */
1990a6053eceSJunchao Zhang   }
19913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19923ab56b82SJunchao Zhang }
19933ab56b82SJunchao Zhang 
1994d2a308c1SPierre Jolivet static PetscErrorCode MatFactorNumeric_MUMPS(Mat F, Mat A, PETSC_UNUSED const MatFactorInfo *info)
1995d71ae5a4SJacob Faibussowitsch {
199657508eceSPierre Jolivet   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
1997ace3abfcSBarry Smith   PetscBool  isMPIAIJ;
1998397b6df1SKris Buschelman 
1999397b6df1SKris Buschelman   PetscFunctionBegin;
2000dbf6bb8dSprj-   if (mumps->id.INFOG(1) < 0 && !(mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0)) {
200148a46eb9SPierre 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)));
20029566063dSJacob 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)));
20033ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
20042aca8efcSHong Zhang   }
20056baea169SHong Zhang 
20069566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, mumps));
20079566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_REUSE_MATRIX, mumps));
2008397b6df1SKris Buschelman 
2009397b6df1SKris Buschelman   /* numerical factorization phase */
2010a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
20114e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
2012ad540459SPierre Jolivet     if (!mumps->myid) mumps->id.a = (MumpsScalar *)mumps->val;
2013397b6df1SKris Buschelman   } else {
2014940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar *)mumps->val;
2015397b6df1SKris Buschelman   }
20163ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
2017a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
20189261f6e4SBarry Smith     PetscCheck(!A->erroriffailure, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in numerical factorization: INFOG(1)=%d, INFO(2)=%d " MUMPS_MANUALS, mumps->id.INFOG(1), mumps->id.INFO(2));
20199261f6e4SBarry Smith     if (mumps->id.INFOG(1) == -10) {
20209261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in numerical factorization: matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2021603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2022c0d63f2fSHong Zhang     } else if (mumps->id.INFOG(1) == -13) {
20239261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in numerical factorization: INFOG(1)=%d, cannot allocate required memory %d megabytes\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2024603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OUTMEMORY;
2025c0d63f2fSHong Zhang     } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10)) {
2026bdcd51b8SPierre Jolivet       PetscCall(PetscInfo(F, "MUMPS error in numerical factorization: INFOG(1)=%d, INFO(2)=%d, problem with work array\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2027603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OUTMEMORY;
20282aca8efcSHong Zhang     } else {
20299261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in numerical factorization: INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2030603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OTHER;
2031151787a6SHong Zhang     }
20322aca8efcSHong Zhang   }
20339261f6e4SBarry Smith   PetscCheck(mumps->myid || mumps->id.ICNTL(16) <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in numerical factorization: ICNTL(16)=%d " MUMPS_MANUALS, mumps->id.INFOG(16));
2034397b6df1SKris Buschelman 
2035b3cb21ddSStefano Zampini   F->assembled = PETSC_TRUE;
2036d47f36abSHong Zhang 
2037b3cb21ddSStefano Zampini   if (F->schur) { /* reset Schur status to unfactored */
20383cb7dd0eSStefano Zampini #if defined(PETSC_HAVE_CUDA)
2039c70f7ee4SJunchao Zhang     F->schur->offloadmask = PETSC_OFFLOAD_CPU;
20403cb7dd0eSStefano Zampini #endif
2041b3cb21ddSStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
2042b3cb21ddSStefano Zampini       mumps->id.ICNTL(19) = 2;
20439566063dSJacob Faibussowitsch       PetscCall(MatTranspose(F->schur, MAT_INPLACE_MATRIX, &F->schur));
2044b3cb21ddSStefano Zampini     }
20459566063dSJacob Faibussowitsch     PetscCall(MatFactorRestoreSchurComplement(F, NULL, MAT_FACTOR_SCHUR_UNFACTORED));
2046b3cb21ddSStefano Zampini   }
204767877ebaSShri Abhyankar 
2048066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
2049066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
2050066565c5SStefano Zampini 
20513ab56b82SJunchao Zhang   if (!mumps->is_omp_master) mumps->id.INFO(23) = 0;
20522d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
205367877ebaSShri Abhyankar     PetscInt     lsol_loc;
205467877ebaSShri Abhyankar     PetscScalar *sol_loc;
20552205254eSKarl Rupp 
20569566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ));
2057c2093ab7SHong Zhang 
2058c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
2059c2093ab7SHong Zhang     if (mumps->x_seq) {
20609566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&mumps->scat_sol));
20619566063dSJacob Faibussowitsch       PetscCall(PetscFree2(mumps->id.sol_loc, mumps->id.isol_loc));
20629566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mumps->x_seq));
2063c2093ab7SHong Zhang     }
2064a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
20659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lsol_loc, &sol_loc, lsol_loc, &mumps->id.isol_loc));
20666497c311SBarry Smith     mumps->id.lsol_loc = (PetscMUMPSInt)lsol_loc;
2067940cd9d6SSatish Balay     mumps->id.sol_loc  = (MumpsScalar *)sol_loc;
20689566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, lsol_loc, sol_loc, &mumps->x_seq));
206967877ebaSShri Abhyankar   }
2070d2a308c1SPierre Jolivet   PetscCall(PetscLogFlops((double)mumps->id.RINFO(2)));
20713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2072397b6df1SKris Buschelman }
2073397b6df1SKris Buschelman 
20749a2535b5SHong Zhang /* Sets MUMPS options from the options database */
207566976f2fSJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_MUMPS(Mat F, Mat A)
2076d71ae5a4SJacob Faibussowitsch {
2077e69c285eSBarry Smith   Mat_MUMPS    *mumps = (Mat_MUMPS *)F->data;
2078413bcc21SPierre Jolivet   PetscMUMPSInt icntl = 0, size, *listvar_schur;
207945e3843bSPierre Jolivet   PetscInt      info[80], i, ninfo = 80, rbs, cbs;
2080413bcc21SPierre Jolivet   PetscBool     flg = PETSC_FALSE, schur = (PetscBool)(mumps->id.ICNTL(26) == -1);
2081413bcc21SPierre Jolivet   MumpsScalar  *arr;
2082dcd589f8SShri Abhyankar 
2083dcd589f8SShri Abhyankar   PetscFunctionBegin;
208426cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "MUMPS Options", "Mat");
2085413bcc21SPierre Jolivet   if (mumps->id.job == JOB_NULL) { /* MatSetFromOptions_MUMPS() has never been called before */
2086413bcc21SPierre Jolivet     PetscInt nthreads   = 0;
2087413bcc21SPierre Jolivet     PetscInt nCNTL_pre  = mumps->CNTL_pre ? mumps->CNTL_pre[0] : 0;
2088413bcc21SPierre Jolivet     PetscInt nICNTL_pre = mumps->ICNTL_pre ? mumps->ICNTL_pre[0] : 0;
2089413bcc21SPierre Jolivet 
2090413bcc21SPierre Jolivet     mumps->petsc_comm = PetscObjectComm((PetscObject)A);
2091413bcc21SPierre Jolivet     PetscCallMPI(MPI_Comm_size(mumps->petsc_comm, &mumps->petsc_size));
2092413bcc21SPierre Jolivet     PetscCallMPI(MPI_Comm_rank(mumps->petsc_comm, &mumps->myid)); /* "if (!myid)" still works even if mumps_comm is different */
2093413bcc21SPierre Jolivet 
2094413bcc21SPierre Jolivet     PetscCall(PetscOptionsName("-mat_mumps_use_omp_threads", "Convert MPI processes into OpenMP threads", "None", &mumps->use_petsc_omp_support));
2095413bcc21SPierre Jolivet     if (mumps->use_petsc_omp_support) nthreads = -1; /* -1 will let PetscOmpCtrlCreate() guess a proper value when user did not supply one */
2096413bcc21SPierre Jolivet     /* do not use PetscOptionsInt() so that the option -mat_mumps_use_omp_threads is not displayed twice in the help */
2097413bcc21SPierre Jolivet     PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)F)->prefix, "-mat_mumps_use_omp_threads", &nthreads, NULL));
2098413bcc21SPierre Jolivet     if (mumps->use_petsc_omp_support) {
2099413bcc21SPierre 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 : "");
2100413bcc21SPierre Jolivet #if defined(PETSC_HAVE_OPENMP_SUPPORT)
2101413bcc21SPierre Jolivet       PetscCall(PetscOmpCtrlCreate(mumps->petsc_comm, nthreads, &mumps->omp_ctrl));
2102413bcc21SPierre Jolivet       PetscCall(PetscOmpCtrlGetOmpComms(mumps->omp_ctrl, &mumps->omp_comm, &mumps->mumps_comm, &mumps->is_omp_master));
2103ea17275aSJose E. Roman #else
2104ea17275aSJose E. Roman       SETERRQ(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",
2105ea17275aSJose E. Roman               ((PetscObject)F)->prefix ? ((PetscObject)F)->prefix : "");
2106413bcc21SPierre Jolivet #endif
2107413bcc21SPierre Jolivet     } else {
2108413bcc21SPierre Jolivet       mumps->omp_comm      = PETSC_COMM_SELF;
2109413bcc21SPierre Jolivet       mumps->mumps_comm    = mumps->petsc_comm;
2110413bcc21SPierre Jolivet       mumps->is_omp_master = PETSC_TRUE;
2111413bcc21SPierre Jolivet     }
2112413bcc21SPierre Jolivet     PetscCallMPI(MPI_Comm_size(mumps->omp_comm, &mumps->omp_comm_size));
2113413bcc21SPierre Jolivet     mumps->reqs = NULL;
2114413bcc21SPierre Jolivet     mumps->tag  = 0;
2115413bcc21SPierre Jolivet 
2116413bcc21SPierre Jolivet     if (mumps->mumps_comm != MPI_COMM_NULL) {
2117413bcc21SPierre Jolivet       if (PetscDefined(HAVE_OPENMP_SUPPORT) && mumps->use_petsc_omp_support) {
2118413bcc21SPierre Jolivet         /* It looks like MUMPS does not dup the input comm. Dup a new comm for MUMPS to avoid any tag mismatches. */
2119413bcc21SPierre Jolivet         MPI_Comm comm;
2120413bcc21SPierre Jolivet         PetscCallMPI(MPI_Comm_dup(mumps->mumps_comm, &comm));
2121413bcc21SPierre Jolivet         mumps->mumps_comm = comm;
2122413bcc21SPierre Jolivet       } else PetscCall(PetscCommGetComm(mumps->petsc_comm, &mumps->mumps_comm));
2123413bcc21SPierre Jolivet     }
2124413bcc21SPierre Jolivet 
2125413bcc21SPierre Jolivet     mumps->id.comm_fortran = MPI_Comm_c2f(mumps->mumps_comm);
2126413bcc21SPierre Jolivet     mumps->id.job          = JOB_INIT;
2127413bcc21SPierre Jolivet     mumps->id.par          = 1; /* host participates factorizaton and solve */
2128413bcc21SPierre Jolivet     mumps->id.sym          = mumps->sym;
2129413bcc21SPierre Jolivet 
2130413bcc21SPierre Jolivet     size          = mumps->id.size_schur;
2131413bcc21SPierre Jolivet     arr           = mumps->id.schur;
2132413bcc21SPierre Jolivet     listvar_schur = mumps->id.listvar_schur;
2133413bcc21SPierre Jolivet     PetscMUMPS_c(mumps);
21349261f6e4SBarry Smith     PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
213551ad14ebSPierre Jolivet 
213651ad14ebSPierre Jolivet     /* set PETSc-MUMPS default options - override MUMPS default */
213751ad14ebSPierre Jolivet     mumps->id.ICNTL(3) = 0;
213851ad14ebSPierre Jolivet     mumps->id.ICNTL(4) = 0;
213951ad14ebSPierre Jolivet     if (mumps->petsc_size == 1) {
214051ad14ebSPierre Jolivet       mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */
214151ad14ebSPierre Jolivet       mumps->id.ICNTL(7)  = 7; /* automatic choice of ordering done by the package */
214251ad14ebSPierre Jolivet     } else {
214351ad14ebSPierre Jolivet       mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */
214451ad14ebSPierre Jolivet       mumps->id.ICNTL(21) = 1; /* distributed solution */
214551ad14ebSPierre Jolivet     }
214651ad14ebSPierre Jolivet 
2147413bcc21SPierre Jolivet     /* restore cached ICNTL and CNTL values */
2148413bcc21SPierre Jolivet     for (icntl = 0; icntl < nICNTL_pre; ++icntl) mumps->id.ICNTL(mumps->ICNTL_pre[1 + 2 * icntl]) = mumps->ICNTL_pre[2 + 2 * icntl];
2149413bcc21SPierre Jolivet     for (icntl = 0; icntl < nCNTL_pre; ++icntl) mumps->id.CNTL((PetscInt)mumps->CNTL_pre[1 + 2 * icntl]) = mumps->CNTL_pre[2 + 2 * icntl];
2150413bcc21SPierre Jolivet     PetscCall(PetscFree(mumps->ICNTL_pre));
2151413bcc21SPierre Jolivet     PetscCall(PetscFree(mumps->CNTL_pre));
2152413bcc21SPierre Jolivet 
2153413bcc21SPierre Jolivet     if (schur) {
2154413bcc21SPierre Jolivet       mumps->id.size_schur    = size;
2155413bcc21SPierre Jolivet       mumps->id.schur_lld     = size;
2156413bcc21SPierre Jolivet       mumps->id.schur         = arr;
2157413bcc21SPierre Jolivet       mumps->id.listvar_schur = listvar_schur;
2158413bcc21SPierre Jolivet       if (mumps->petsc_size > 1) {
2159413bcc21SPierre Jolivet         PetscBool gs; /* gs is false if any rank other than root has non-empty IS */
2160413bcc21SPierre Jolivet 
2161413bcc21SPierre Jolivet         mumps->id.ICNTL(19) = 1;                                                                            /* MUMPS returns Schur centralized on the host */
2162413bcc21SPierre 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 */
2163462c564dSBarry Smith         PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &gs, 1, MPIU_BOOL, MPI_LAND, mumps->petsc_comm));
2164413bcc21SPierre Jolivet         PetscCheck(gs, PETSC_COMM_SELF, PETSC_ERR_SUP, "MUMPS distributed parallel Schur complements not yet supported from PETSc");
2165413bcc21SPierre Jolivet       } else {
2166413bcc21SPierre Jolivet         if (F->factortype == MAT_FACTOR_LU) {
2167413bcc21SPierre Jolivet           mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
2168413bcc21SPierre Jolivet         } else {
2169413bcc21SPierre Jolivet           mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
2170413bcc21SPierre Jolivet         }
2171413bcc21SPierre Jolivet       }
2172413bcc21SPierre Jolivet       mumps->id.ICNTL(26) = -1;
2173413bcc21SPierre Jolivet     }
2174413bcc21SPierre Jolivet 
2175413bcc21SPierre Jolivet     /* copy MUMPS default control values from master to slaves. Although slaves do not call MUMPS, they may access these values in code.
2176413bcc21SPierre Jolivet        For example, ICNTL(9) is initialized to 1 by MUMPS and slaves check ICNTL(9) in MatSolve_MUMPS.
2177413bcc21SPierre Jolivet      */
2178413bcc21SPierre Jolivet     PetscCallMPI(MPI_Bcast(mumps->id.icntl, 40, MPI_INT, 0, mumps->omp_comm));
2179413bcc21SPierre Jolivet     PetscCallMPI(MPI_Bcast(mumps->id.cntl, 15, MPIU_REAL, 0, mumps->omp_comm));
2180413bcc21SPierre Jolivet 
2181413bcc21SPierre Jolivet     mumps->scat_rhs = NULL;
2182413bcc21SPierre Jolivet     mumps->scat_sol = NULL;
2183413bcc21SPierre Jolivet   }
21849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_1", "ICNTL(1): output stream for error messages", "None", mumps->id.ICNTL(1), &icntl, &flg));
21859a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
21869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_2", "ICNTL(2): output stream for diagnostic printing, statistics, and warning", "None", mumps->id.ICNTL(2), &icntl, &flg));
21879a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
21889566063dSJacob 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));
21899a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
2190dcd589f8SShri Abhyankar 
21919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_4", "ICNTL(4): level of printing (0 to 4)", "None", mumps->id.ICNTL(4), &icntl, &flg));
21929a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
21939a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
21949a2535b5SHong Zhang 
21959566063dSJacob 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));
21969a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
21979a2535b5SHong Zhang 
21989566063dSJacob 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));
2199dcd589f8SShri Abhyankar   if (flg) {
2200aed4548fSBarry 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");
2201b53c1a7fSBarry Smith     mumps->id.ICNTL(7) = icntl;
2202dcd589f8SShri Abhyankar   }
2203e0b74bf9SHong Zhang 
22049566063dSJacob 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));
22059566063dSJacob 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() */
22069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_10", "ICNTL(10): max num of refinements", "None", mumps->id.ICNTL(10), &mumps->id.ICNTL(10), NULL));
22079566063dSJacob 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));
22089566063dSJacob 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));
22099566063dSJacob 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));
22109566063dSJacob 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));
221145e3843bSPierre Jolivet   PetscCall(MatGetBlockSizes(A, &rbs, &cbs));
22126497c311SBarry Smith   if (rbs == cbs && rbs > 1) mumps->id.ICNTL(15) = (PetscMUMPSInt)-rbs;
221345e3843bSPierre 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));
221445e3843bSPierre Jolivet   if (flg) {
221545e3843bSPierre Jolivet     PetscCheck(mumps->id.ICNTL(15) <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Positive -mat_mumps_icntl_15 not handled");
221645e3843bSPierre 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");
221745e3843bSPierre Jolivet   }
22189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_19", "ICNTL(19): computes the Schur complement", "None", mumps->id.ICNTL(19), &mumps->id.ICNTL(19), NULL));
221959ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
22209566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F->schur));
22219566063dSJacob Faibussowitsch     PetscCall(MatMumpsResetSchur_Private(mumps));
222259ac8732SStefano Zampini   }
222325aac85cSJunchao Zhang 
222443f3b051SJunchao Zhang   /* Two MPICH Fortran MPI_IN_PLACE binding bugs prevented the use of 'mpich + mumps'. One happened with "mpi4py + mpich + mumps",
222543f3b051SJunchao Zhang      and was reported by Firedrake. See https://bitbucket.org/mpi4py/mpi4py/issues/162/mpi4py-initialization-breaks-fortran
222625aac85cSJunchao Zhang      and a petsc-maint mailing list thread with subject 'MUMPS segfaults in parallel because of ...'
222743f3b051SJunchao Zhang      This bug was fixed by https://github.com/pmodels/mpich/pull/4149. But the fix brought a new bug,
222843f3b051SJunchao Zhang      see https://github.com/pmodels/mpich/issues/5589. This bug was fixed by https://github.com/pmodels/mpich/pull/5590.
222941caa250SJunchao Zhang      In short, we could not use distributed RHS until with MPICH v4.0b1 or we enabled a workaround in mumps-5.6.2+
223025aac85cSJunchao Zhang    */
223141caa250SJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_GE(5, 6, 2) && defined(PETSC_HAVE_MUMPS_AVOID_MPI_IN_PLACE)
223241caa250SJunchao Zhang   mumps->ICNTL20 = 10;
2233e3c15826SJunchao Zhang #elif PETSC_PKG_MUMPS_VERSION_LT(5, 3, 0) || (defined(PETSC_HAVE_MPICH) && (MPICH_NUMVERSION < 40000101))
223425aac85cSJunchao Zhang   mumps->ICNTL20 = 0; /* Centralized dense RHS*/
223543f3b051SJunchao Zhang #else
223643f3b051SJunchao Zhang   mumps->ICNTL20 = 10; /* Distributed dense RHS*/
223725aac85cSJunchao Zhang #endif
22389566063dSJacob 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));
2239aed4548fSBarry 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);
224025aac85cSJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_LT(5, 3, 0)
2241aed4548fSBarry Smith   PetscCheck(!flg || mumps->ICNTL20 != 10, PETSC_COMM_SELF, PETSC_ERR_SUP, "ICNTL(20)=10 is not supported before MUMPS-5.3.0");
224225aac85cSJunchao Zhang #endif
22439566063dSJacob 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 */
22449a2535b5SHong Zhang 
22459566063dSJacob 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));
22469566063dSJacob 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));
22479566063dSJacob 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));
22489371c9d4SSatish Balay   if (mumps->id.ICNTL(24)) { mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ }
2249d7ebd59bSHong Zhang 
22509566063dSJacob 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));
22519566063dSJacob 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));
22529566063dSJacob 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));
2253fa6fd9d0SPierre Jolivet   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_28", "ICNTL(28): use 1 for sequential analysis and ICNTL(7) ordering, or 2 for parallel analysis and ICNTL(29) ordering", "None", mumps->id.ICNTL(28), &mumps->id.ICNTL(28), NULL));
22549566063dSJacob 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));
22559566063dSJacob 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 */
22569566063dSJacob 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));
2257145b44c9SPierre Jolivet   /* PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elimination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL));  -- not supported by PETSc API */
22589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_33", "ICNTL(33): compute determinant", "None", mumps->id.ICNTL(33), &mumps->id.ICNTL(33), NULL));
22599566063dSJacob 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));
22609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_36", "ICNTL(36): choice of BLR factorization variant", "None", mumps->id.ICNTL(36), &mumps->id.ICNTL(36), NULL));
2261*50ea2040Saszaboa   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_37", "ICNTL(37): compression of the contribution blocks (CB)", "None", mumps->id.ICNTL(37), &mumps->id.ICNTL(37), NULL));
22629566063dSJacob 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));
2263c92b4f89SPierre Jolivet   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_48", "ICNTL(48): multithreading with tree parallelism", "None", mumps->id.ICNTL(48), &mumps->id.ICNTL(48), NULL));
2264146931dbSPierre Jolivet   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_58", "ICNTL(58): defines options for symbolic factorization", "None", mumps->id.ICNTL(58), &mumps->id.ICNTL(58), NULL));
2265dcd589f8SShri Abhyankar 
22669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_1", "CNTL(1): relative pivoting threshold", "None", mumps->id.CNTL(1), &mumps->id.CNTL(1), NULL));
22679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_2", "CNTL(2): stopping criterion of refinement", "None", mumps->id.CNTL(2), &mumps->id.CNTL(2), NULL));
22689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_3", "CNTL(3): absolute pivoting threshold", "None", mumps->id.CNTL(3), &mumps->id.CNTL(3), NULL));
22699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_4", "CNTL(4): value for static pivoting", "None", mumps->id.CNTL(4), &mumps->id.CNTL(4), NULL));
22709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_5", "CNTL(5): fixation for null pivots", "None", mumps->id.CNTL(5), &mumps->id.CNTL(5), NULL));
22719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_7", "CNTL(7): dropping parameter used during BLR", "None", mumps->id.CNTL(7), &mumps->id.CNTL(7), NULL));
2272e5bb22a1SHong Zhang 
22739566063dSJacob 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));
2274b34f08ffSHong Zhang 
22759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-mat_mumps_view_info", "request INFO local to each processor", "", info, &ninfo, NULL));
2276b34f08ffSHong Zhang   if (ninfo) {
227708401ef6SPierre Jolivet     PetscCheck(ninfo <= 80, PETSC_COMM_SELF, PETSC_ERR_USER, "number of INFO %" PetscInt_FMT " must <= 80", ninfo);
22789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ninfo, &mumps->info));
2279b34f08ffSHong Zhang     mumps->ninfo = ninfo;
2280b34f08ffSHong Zhang     for (i = 0; i < ninfo; i++) {
2281aed4548fSBarry 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);
2282f7d195e4SLawrence Mitchell       mumps->info[i] = info[i];
2283b34f08ffSHong Zhang     }
2284b34f08ffSHong Zhang   }
2285d0609cedSBarry Smith   PetscOptionsEnd();
22863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2287dcd589f8SShri Abhyankar }
2288dcd589f8SShri Abhyankar 
2289d2a308c1SPierre Jolivet static PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F, Mat A, PETSC_UNUSED const MatFactorInfo *info, Mat_MUMPS *mumps)
2290d71ae5a4SJacob Faibussowitsch {
22915cd7cf9dSHong Zhang   PetscFunctionBegin;
22925cd7cf9dSHong Zhang   if (mumps->id.INFOG(1) < 0) {
22939261f6e4SBarry Smith     PetscCheck(!A->erroriffailure, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in analysis: INFOG(1)=%d " MUMPS_MANUALS, mumps->id.INFOG(1));
22945cd7cf9dSHong Zhang     if (mumps->id.INFOG(1) == -6) {
22959261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in analysis: matrix is singular, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2296603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
22975cd7cf9dSHong Zhang     } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
22989261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in analysis: problem with work array, INFOG(1)=%d, INFO(2)=%d\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2299603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OUTMEMORY;
2300dbf6bb8dSprj-     } else if (mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0) {
23019261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in analysis: empty matrix\n"));
23025cd7cf9dSHong Zhang     } else {
23039261f6e4SBarry Smith       PetscCall(PetscInfo(F, "MUMPS error in analysis: INFOG(1)=%d, INFO(2)=%d " MUMPS_MANUALS "\n", mumps->id.INFOG(1), mumps->id.INFO(2)));
2304603e8f96SBarry Smith       F->factorerrortype = MAT_FACTOR_OTHER;
23055cd7cf9dSHong Zhang     }
23065cd7cf9dSHong Zhang   }
230772b150d8SStefano Zampini   if (!mumps->id.n) F->factorerrortype = MAT_FACTOR_NOERROR;
23083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23095cd7cf9dSHong Zhang }
23105cd7cf9dSHong Zhang 
2311d2a308c1SPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F, Mat A, IS r, PETSC_UNUSED IS c, const MatFactorInfo *info)
2312d71ae5a4SJacob Faibussowitsch {
2313e69c285eSBarry Smith   Mat_MUMPS     *mumps = (Mat_MUMPS *)F->data;
231467877ebaSShri Abhyankar   Vec            b;
231567877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
2316397b6df1SKris Buschelman 
2317397b6df1SKris Buschelman   PetscFunctionBegin;
2318d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2319d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
23203ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2321d47f36abSHong Zhang   }
2322dcd589f8SShri Abhyankar 
23239a2535b5SHong Zhang   /* Set MUMPS options from the options database */
232426cc229bSBarry Smith   PetscCall(MatSetFromOptions_MUMPS(F, A));
2325dcd589f8SShri Abhyankar 
23269566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
23279566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX, mumps));
2328dcd589f8SShri Abhyankar 
232967877ebaSShri Abhyankar   /* analysis phase */
2330a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
23316497c311SBarry Smith   PetscCall(PetscMUMPSIntCast(M, &mumps->id.n));
2332a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
233367877ebaSShri Abhyankar   case 0: /* centralized assembled matrix input */
2334a5e57a09SHong Zhang     if (!mumps->myid) {
2335a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2336a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2337a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2338a6053eceSJunchao Zhang       if (mumps->id.ICNTL(6) > 1) mumps->id.a = (MumpsScalar *)mumps->val;
233951ad14ebSPierre Jolivet       if (r && mumps->id.ICNTL(7) == 7) {
23404ac6704cSBarry Smith         mumps->id.ICNTL(7) = 1;
2341a5e57a09SHong Zhang         if (!mumps->myid) {
2342e0b74bf9SHong Zhang           const PetscInt *idx;
2343a6053eceSJunchao Zhang           PetscInt        i;
23442205254eSKarl Rupp 
23459566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(M, &mumps->id.perm_in));
23469566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(r, &idx));
2347f4f49eeaSPierre Jolivet           for (i = 0; i < M; i++) PetscCall(PetscMUMPSIntCast(idx[i] + 1, &mumps->id.perm_in[i])); /* perm_in[]: start from 1, not 0! */
23489566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(r, &idx));
2349e0b74bf9SHong Zhang         }
2350e0b74bf9SHong Zhang       }
235167877ebaSShri Abhyankar     }
235267877ebaSShri Abhyankar     break;
235367877ebaSShri Abhyankar   case 3: /* distributed assembled matrix input (size>1) */
2354a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2355a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2356a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2357a6053eceSJunchao Zhang     if (mumps->id.ICNTL(6) > 1) mumps->id.a_loc = (MumpsScalar *)mumps->val;
235825aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
23599566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
23609566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b, &mumps->scat_rhs, &mumps->b_seq));
23619566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
236225aac85cSJunchao Zhang     }
236367877ebaSShri Abhyankar     break;
236467877ebaSShri Abhyankar   }
23653ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
23669566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F, A, info, mumps));
236767877ebaSShri Abhyankar 
2368719d5645SBarry Smith   F->ops->lufactornumeric   = MatFactorNumeric_MUMPS;
2369dcd589f8SShri Abhyankar   F->ops->solve             = MatSolve_MUMPS;
237051d5961aSHong Zhang   F->ops->solvetranspose    = MatSolveTranspose_MUMPS;
23714e34a73bSHong Zhang   F->ops->matsolve          = MatMatSolve_MUMPS;
2372eb3ef3b2SHong Zhang   F->ops->mattransposesolve = MatMatTransposeSolve_MUMPS;
2373b18964edSHong Zhang   F->ops->matsolvetranspose = MatMatSolveTranspose_MUMPS;
2374d47f36abSHong Zhang 
2375d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
23763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2377b24902e0SBarry Smith }
2378b24902e0SBarry Smith 
2379450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
2380d2a308c1SPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F, Mat A, PETSC_UNUSED IS r, PETSC_UNUSED IS c, const MatFactorInfo *info)
2381d71ae5a4SJacob Faibussowitsch {
2382e69c285eSBarry Smith   Mat_MUMPS     *mumps = (Mat_MUMPS *)F->data;
238367877ebaSShri Abhyankar   Vec            b;
238467877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
2385450b117fSShri Abhyankar 
2386450b117fSShri Abhyankar   PetscFunctionBegin;
2387d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2388338d3105SPierre Jolivet     /* F is assembled by a previous call of MatLUFactorSymbolic_BAIJMUMPS() */
23893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2390d47f36abSHong Zhang   }
2391dcd589f8SShri Abhyankar 
23929a2535b5SHong Zhang   /* Set MUMPS options from the options database */
239326cc229bSBarry Smith   PetscCall(MatSetFromOptions_MUMPS(F, A));
2394dcd589f8SShri Abhyankar 
23959566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
23969566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX, mumps));
239767877ebaSShri Abhyankar 
239867877ebaSShri Abhyankar   /* analysis phase */
2399a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
24006497c311SBarry Smith   PetscCall(PetscMUMPSIntCast(M, &mumps->id.n));
2401a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
240267877ebaSShri Abhyankar   case 0: /* centralized assembled matrix input */
2403a5e57a09SHong Zhang     if (!mumps->myid) {
2404a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2405a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2406a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2407ad540459SPierre Jolivet       if (mumps->id.ICNTL(6) > 1) mumps->id.a = (MumpsScalar *)mumps->val;
240867877ebaSShri Abhyankar     }
240967877ebaSShri Abhyankar     break;
241067877ebaSShri Abhyankar   case 3: /* distributed assembled matrix input (size>1) */
2411a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2412a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2413a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2414ad540459SPierre Jolivet     if (mumps->id.ICNTL(6) > 1) mumps->id.a_loc = (MumpsScalar *)mumps->val;
241525aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
24169566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
24179566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b, &mumps->scat_rhs, &mumps->b_seq));
24189566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
241925aac85cSJunchao Zhang     }
242067877ebaSShri Abhyankar     break;
242167877ebaSShri Abhyankar   }
24223ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
24239566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F, A, info, mumps));
242467877ebaSShri Abhyankar 
2425450b117fSShri Abhyankar   F->ops->lufactornumeric   = MatFactorNumeric_MUMPS;
2426dcd589f8SShri Abhyankar   F->ops->solve             = MatSolve_MUMPS;
242751d5961aSHong Zhang   F->ops->solvetranspose    = MatSolveTranspose_MUMPS;
2428b18964edSHong Zhang   F->ops->matsolvetranspose = MatMatSolveTranspose_MUMPS;
2429d47f36abSHong Zhang 
2430d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
24313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2432450b117fSShri Abhyankar }
2433b24902e0SBarry Smith 
2434141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
2435d2a308c1SPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F, Mat A, PETSC_UNUSED IS r, const MatFactorInfo *info)
2436d71ae5a4SJacob Faibussowitsch {
2437e69c285eSBarry Smith   Mat_MUMPS     *mumps = (Mat_MUMPS *)F->data;
243867877ebaSShri Abhyankar   Vec            b;
243967877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
2440397b6df1SKris Buschelman 
2441397b6df1SKris Buschelman   PetscFunctionBegin;
2442d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2443338d3105SPierre Jolivet     /* F is assembled by a previous call of MatCholeskyFactorSymbolic_MUMPS() */
24443ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2445d47f36abSHong Zhang   }
2446dcd589f8SShri Abhyankar 
24479a2535b5SHong Zhang   /* Set MUMPS options from the options database */
244826cc229bSBarry Smith   PetscCall(MatSetFromOptions_MUMPS(F, A));
2449dcd589f8SShri Abhyankar 
24509566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
24519566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX, mumps));
2452dcd589f8SShri Abhyankar 
245367877ebaSShri Abhyankar   /* analysis phase */
2454a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
24556497c311SBarry Smith   PetscCall(PetscMUMPSIntCast(M, &mumps->id.n));
2456a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
245767877ebaSShri Abhyankar   case 0: /* centralized assembled matrix input */
2458a5e57a09SHong Zhang     if (!mumps->myid) {
2459a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2460a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2461a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2462ad540459SPierre Jolivet       if (mumps->id.ICNTL(6) > 1) mumps->id.a = (MumpsScalar *)mumps->val;
246367877ebaSShri Abhyankar     }
246467877ebaSShri Abhyankar     break;
246567877ebaSShri Abhyankar   case 3: /* distributed assembled matrix input (size>1) */
2466a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2467a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2468a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2469ad540459SPierre Jolivet     if (mumps->id.ICNTL(6) > 1) mumps->id.a_loc = (MumpsScalar *)mumps->val;
247025aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
24719566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
24729566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b, &mumps->scat_rhs, &mumps->b_seq));
24739566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
247425aac85cSJunchao Zhang     }
247567877ebaSShri Abhyankar     break;
247667877ebaSShri Abhyankar   }
24773ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
24789566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F, A, info, mumps));
24795cd7cf9dSHong Zhang 
24802792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
2481dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
248251d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
24834e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
248423a5080aSHong Zhang   F->ops->mattransposesolve     = MatMatTransposeSolve_MUMPS;
2485b18964edSHong Zhang   F->ops->matsolvetranspose     = MatMatSolveTranspose_MUMPS;
24864e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
24870298fd71SBarry Smith   F->ops->getinertia = NULL;
24884e34a73bSHong Zhang #else
24894e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
2490db4efbfdSBarry Smith #endif
2491d47f36abSHong Zhang 
2492d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
24933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2494b24902e0SBarry Smith }
2495b24902e0SBarry Smith 
249666976f2fSJacob Faibussowitsch static PetscErrorCode MatView_MUMPS(Mat A, PetscViewer viewer)
2497d71ae5a4SJacob Faibussowitsch {
249864e6c443SBarry Smith   PetscBool         iascii;
249964e6c443SBarry Smith   PetscViewerFormat format;
2500e69c285eSBarry Smith   Mat_MUMPS        *mumps = (Mat_MUMPS *)A->data;
2501f6c57405SHong Zhang 
2502f6c57405SHong Zhang   PetscFunctionBegin;
250364e6c443SBarry Smith   /* check if matrix is mumps type */
25043ba16761SJacob Faibussowitsch   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(PETSC_SUCCESS);
250564e6c443SBarry Smith 
25069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
250764e6c443SBarry Smith   if (iascii) {
25089566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
25091511cd71SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
25109566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "MUMPS run parameters:\n"));
25111511cd71SPierre Jolivet       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
25129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  SYM (matrix type):                   %d\n", mumps->id.sym));
25139566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  PAR (host participation):            %d\n", mumps->id.par));
25149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(1) (output for error):         %d\n", mumps->id.ICNTL(1)));
25159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(2) (output of diagnostic msg): %d\n", mumps->id.ICNTL(2)));
25169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(3) (output for global info):   %d\n", mumps->id.ICNTL(3)));
25179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(4) (level of printing):        %d\n", mumps->id.ICNTL(4)));
25189566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(5) (input mat struct):         %d\n", mumps->id.ICNTL(5)));
25199566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(6) (matrix prescaling):        %d\n", mumps->id.ICNTL(6)));
25209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(7) (sequential matrix ordering):%d\n", mumps->id.ICNTL(7)));
25219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(8) (scaling strategy):         %d\n", mumps->id.ICNTL(8)));
25229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(10) (max num of refinements):  %d\n", mumps->id.ICNTL(10)));
25239566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(11) (error analysis):          %d\n", mumps->id.ICNTL(11)));
2524a5e57a09SHong Zhang         if (mumps->id.ICNTL(11) > 0) {
2525d2a308c1SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(4) (inf norm of input mat):        %g\n", (double)mumps->id.RINFOG(4)));
2526d2a308c1SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(5) (inf norm of solution):         %g\n", (double)mumps->id.RINFOG(5)));
2527d2a308c1SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(6) (inf norm of residual):         %g\n", (double)mumps->id.RINFOG(6)));
2528d2a308c1SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n", (double)mumps->id.RINFOG(7), (double)mumps->id.RINFOG(8)));
2529d2a308c1SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(9) (error estimate):               %g\n", (double)mumps->id.RINFOG(9)));
2530d2a308c1SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n", (double)mumps->id.RINFOG(10), (double)mumps->id.RINFOG(11)));
2531f6c57405SHong Zhang         }
25329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(12) (efficiency control):                         %d\n", mumps->id.ICNTL(12)));
25339566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(13) (sequential factorization of the root node):  %d\n", mumps->id.ICNTL(13)));
25349566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(14) (percentage of estimated workspace increase): %d\n", mumps->id.ICNTL(14)));
253545e3843bSPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(15) (compression of the input matrix):            %d\n", mumps->id.ICNTL(15)));
2536f6c57405SHong Zhang         /* ICNTL(15-17) not used */
25379566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(18) (input mat struct):                           %d\n", mumps->id.ICNTL(18)));
25389566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(19) (Schur complement info):                      %d\n", mumps->id.ICNTL(19)));
25399566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(20) (RHS sparse pattern):                         %d\n", mumps->id.ICNTL(20)));
25409566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(21) (solution struct):                            %d\n", mumps->id.ICNTL(21)));
25419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(22) (in-core/out-of-core facility):               %d\n", mumps->id.ICNTL(22)));
25429566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(23) (max size of memory can be allocated locally):%d\n", mumps->id.ICNTL(23)));
2543c0165424SHong Zhang 
25449566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(24) (detection of null pivot rows):               %d\n", mumps->id.ICNTL(24)));
25459566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(25) (computation of a null space basis):          %d\n", mumps->id.ICNTL(25)));
25469566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(26) (Schur options for RHS or solution):          %d\n", mumps->id.ICNTL(26)));
25479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(27) (blocking size for multiple RHS):             %d\n", mumps->id.ICNTL(27)));
25489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(28) (use parallel or sequential ordering):        %d\n", mumps->id.ICNTL(28)));
25499566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(29) (parallel ordering):                          %d\n", mumps->id.ICNTL(29)));
255042179a6aSHong Zhang 
25519566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(30) (user-specified set of entries in inv(A)):    %d\n", mumps->id.ICNTL(30)));
25529566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(31) (factors is discarded in the solve phase):    %d\n", mumps->id.ICNTL(31)));
25539566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(33) (compute determinant):                        %d\n", mumps->id.ICNTL(33)));
25549566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(35) (activate BLR based factorization):           %d\n", mumps->id.ICNTL(35)));
25559566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(36) (choice of BLR factorization variant):        %d\n", mumps->id.ICNTL(36)));
2556*50ea2040Saszaboa         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(37) (compression of the contribution blocks):     %d\n", mumps->id.ICNTL(37)));
25579566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(38) (estimated compression rate of LU factors):   %d\n", mumps->id.ICNTL(38)));
2558c92b4f89SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(48) (multithreading with tree parallelism):       %d\n", mumps->id.ICNTL(48)));
2559146931dbSPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  ICNTL(58) (options for symbolic factorization):         %d\n", mumps->id.ICNTL(58)));
2560f6c57405SHong Zhang 
2561d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(1) (relative pivoting threshold):      %g\n", (double)mumps->id.CNTL(1)));
2562d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(2) (stopping criterion of refinement): %g\n", (double)mumps->id.CNTL(2)));
2563d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(3) (absolute pivoting threshold):      %g\n", (double)mumps->id.CNTL(3)));
2564d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(4) (value of static pivoting):         %g\n", (double)mumps->id.CNTL(4)));
2565d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(5) (fixation for null pivots):         %g\n", (double)mumps->id.CNTL(5)));
2566d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  CNTL(7) (dropping parameter for BLR):       %g\n", (double)mumps->id.CNTL(7)));
2567f6c57405SHong Zhang 
2568a5b23f4aSJose E. Roman         /* information local to each processor */
25699566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis):\n"));
25709566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2571d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %g\n", mumps->myid, (double)mumps->id.RINFO(1)));
25729566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
25739566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization):\n"));
2574d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %g\n", mumps->myid, (double)mumps->id.RINFO(2)));
25759566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
25769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization):\n"));
2577d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %g\n", mumps->myid, (double)mumps->id.RINFO(3)));
25789566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2579f6c57405SHong Zhang 
25809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization):\n"));
25819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(15)));
25829566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2583f6c57405SHong Zhang 
25849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization):\n"));
25859566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(16)));
25869566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2587f6c57405SHong Zhang 
25889566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization):\n"));
25899566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(23)));
25909566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2591b34f08ffSHong Zhang 
2592a0e18203SThibaut Appel         if (mumps->ninfo && mumps->ninfo <= 80) {
2593b34f08ffSHong Zhang           PetscInt i;
2594b34f08ffSHong Zhang           for (i = 0; i < mumps->ninfo; i++) {
25959566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(%" PetscInt_FMT "):\n", mumps->info[i]));
25969566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "    [%d] %d\n", mumps->myid, mumps->id.INFO(mumps->info[i])));
25979566063dSJacob Faibussowitsch             PetscCall(PetscViewerFlush(viewer));
2598b34f08ffSHong Zhang           }
2599b34f08ffSHong Zhang         }
26009566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
26011511cd71SPierre Jolivet       } else PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all processes\n", ((PetscObject)A)->prefix ? ((PetscObject)A)->prefix : ""));
2602f6c57405SHong Zhang 
26031511cd71SPierre Jolivet       if (mumps->myid == 0) { /* information from the host */
2604d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFOG(1) (global estimated flops for the elimination after analysis): %g\n", (double)mumps->id.RINFOG(1)));
2605d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFOG(2) (global estimated flops for the assembly after factorization): %g\n", (double)mumps->id.RINFOG(2)));
2606d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFOG(3) (global estimated flops for the elimination after factorization): %g\n", (double)mumps->id.RINFOG(3)));
2607d2a308c1SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n", (double)mumps->id.RINFOG(12), (double)mumps->id.RINFOG(13), mumps->id.INFOG(34)));
2608f6c57405SHong Zhang 
26099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d\n", mumps->id.INFOG(3)));
26109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d\n", mumps->id.INFOG(4)));
26119566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(5) (estimated maximum front size in the complete tree): %d\n", mumps->id.INFOG(5)));
26129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(6) (number of nodes in the complete tree): %d\n", mumps->id.INFOG(6)));
26139566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(7) (ordering option effectively used after analysis): %d\n", mumps->id.INFOG(7)));
26149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d\n", mumps->id.INFOG(8)));
26159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d\n", mumps->id.INFOG(9)));
26169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(10) (total integer space store the matrix factors after factorization): %d\n", mumps->id.INFOG(10)));
26179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(11) (order of largest frontal matrix after factorization): %d\n", mumps->id.INFOG(11)));
26189566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(12) (number of off-diagonal pivots): %d\n", mumps->id.INFOG(12)));
26199566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(13) (number of delayed pivots after factorization): %d\n", mumps->id.INFOG(13)));
26209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(14) (number of memory compress after factorization): %d\n", mumps->id.INFOG(14)));
26219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(15) (number of steps of iterative refinement after solution): %d\n", mumps->id.INFOG(15)));
26229566063dSJacob 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)));
26239566063dSJacob 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)));
26249566063dSJacob 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)));
26259566063dSJacob 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)));
26269566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(20) (estimated number of entries in the factors): %d\n", mumps->id.INFOG(20)));
26279566063dSJacob 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)));
26289566063dSJacob 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)));
26299566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d\n", mumps->id.INFOG(23)));
26309566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d\n", mumps->id.INFOG(24)));
26319566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d\n", mumps->id.INFOG(25)));
26329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(28) (after factorization: number of null pivots encountered): %d\n", mumps->id.INFOG(28)));
26339566063dSJacob 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)));
26349566063dSJacob 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)));
26359566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(32) (after analysis: type of analysis done): %d\n", mumps->id.INFOG(32)));
26369566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(33) (value used for ICNTL(8)): %d\n", mumps->id.INFOG(33)));
26379566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  INFOG(34) (exponent of the determinant if determinant is requested): %d\n", mumps->id.INFOG(34)));
26389566063dSJacob 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)));
26399566063dSJacob 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)));
26409566063dSJacob 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)));
26419566063dSJacob 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)));
26429566063dSJacob 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)));
2643f6c57405SHong Zhang       }
2644f6c57405SHong Zhang     }
2645cb828f0fSHong Zhang   }
26463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2647f6c57405SHong Zhang }
2648f6c57405SHong Zhang 
2649d2a308c1SPierre Jolivet static PetscErrorCode MatGetInfo_MUMPS(Mat A, PETSC_UNUSED MatInfoType flag, MatInfo *info)
2650d71ae5a4SJacob Faibussowitsch {
2651e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)A->data;
265235bd34faSBarry Smith 
265335bd34faSBarry Smith   PetscFunctionBegin;
265435bd34faSBarry Smith   info->block_size        = 1.0;
265564412097SPierre Jolivet   info->nz_allocated      = mumps->id.INFOG(20) >= 0 ? mumps->id.INFOG(20) : -1000000 * mumps->id.INFOG(20);
265664412097SPierre Jolivet   info->nz_used           = mumps->id.INFOG(20) >= 0 ? mumps->id.INFOG(20) : -1000000 * mumps->id.INFOG(20);
265735bd34faSBarry Smith   info->nz_unneeded       = 0.0;
265835bd34faSBarry Smith   info->assemblies        = 0.0;
265935bd34faSBarry Smith   info->mallocs           = 0.0;
266035bd34faSBarry Smith   info->memory            = 0.0;
266135bd34faSBarry Smith   info->fill_ratio_given  = 0;
266235bd34faSBarry Smith   info->fill_ratio_needed = 0;
266335bd34faSBarry Smith   info->factor_mallocs    = 0;
26643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266535bd34faSBarry Smith }
266635bd34faSBarry Smith 
266766976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
2668d71ae5a4SJacob Faibussowitsch {
2669e69c285eSBarry Smith   Mat_MUMPS         *mumps = (Mat_MUMPS *)F->data;
2670a3d589ffSStefano Zampini   const PetscScalar *arr;
26718e7ba810SStefano Zampini   const PetscInt    *idxs;
26728e7ba810SStefano Zampini   PetscInt           size, i;
26736444a565SStefano Zampini 
26746444a565SStefano Zampini   PetscFunctionBegin;
26759566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &size));
2676b3cb21ddSStefano Zampini   /* Schur complement matrix */
26779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F->schur));
26789566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size, size, NULL, &F->schur));
26799566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(F->schur, &arr));
2680a3d589ffSStefano Zampini   mumps->id.schur = (MumpsScalar *)arr;
26816497c311SBarry Smith   PetscCall(PetscMUMPSIntCast(size, &mumps->id.size_schur));
26826497c311SBarry Smith   PetscCall(PetscMUMPSIntCast(size, &mumps->id.schur_lld));
26839566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(F->schur, &arr));
268448a46eb9SPierre Jolivet   if (mumps->sym == 1) PetscCall(MatSetOption(F->schur, MAT_SPD, PETSC_TRUE));
2685b3cb21ddSStefano Zampini 
2686b3cb21ddSStefano Zampini   /* MUMPS expects Fortran style indices */
26879566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.listvar_schur));
26889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &mumps->id.listvar_schur));
26899566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxs));
2690f4f49eeaSPierre Jolivet   for (i = 0; i < size; i++) PetscCall(PetscMUMPSIntCast(idxs[i] + 1, &mumps->id.listvar_schur[i]));
26919566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxs));
269259ac8732SStefano Zampini   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
2693b5fa320bSStefano Zampini   mumps->id.ICNTL(26) = -1;
26943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26956444a565SStefano Zampini }
269659ac8732SStefano Zampini 
269766976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F, Mat *S)
2698d71ae5a4SJacob Faibussowitsch {
26996444a565SStefano Zampini   Mat          St;
2700e69c285eSBarry Smith   Mat_MUMPS   *mumps = (Mat_MUMPS *)F->data;
27016444a565SStefano Zampini   PetscScalar *array;
27026444a565SStefano Zampini 
27036444a565SStefano Zampini   PetscFunctionBegin;
27049261f6e4SBarry Smith   PetscCheck(mumps->id.ICNTL(19), PetscObjectComm((PetscObject)F), PETSC_ERR_ORDER, "Schur complement mode not selected! Call MatFactorSetSchurIS() to enable it");
27059566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &St));
27069566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(St, PETSC_DECIDE, PETSC_DECIDE, mumps->id.size_schur, mumps->id.size_schur));
27079566063dSJacob Faibussowitsch   PetscCall(MatSetType(St, MATDENSE));
27089566063dSJacob Faibussowitsch   PetscCall(MatSetUp(St));
27099566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(St, &array));
271059ac8732SStefano Zampini   if (!mumps->sym) {                /* MUMPS always return a full matrix */
27116444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
27126444a565SStefano Zampini       PetscInt i, j, N = mumps->id.size_schur;
27136444a565SStefano Zampini       for (i = 0; i < N; i++) {
27146444a565SStefano Zampini         for (j = 0; j < N; j++) {
27156444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
27166444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j];
27176444a565SStefano Zampini #else
27183e5b40d0SPierre Jolivet           PetscScalar val = mumps->id.schur[i * N + j].r + PETSC_i * mumps->id.schur[i * N + j].i;
27196444a565SStefano Zampini #endif
27206444a565SStefano Zampini           array[j * N + i] = val;
27216444a565SStefano Zampini         }
27226444a565SStefano Zampini       }
27236444a565SStefano Zampini     } else { /* stored by columns */
27249566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array, mumps->id.schur, mumps->id.size_schur * mumps->id.size_schur));
27256444a565SStefano Zampini     }
27266444a565SStefano Zampini   } else {                          /* either full or lower-triangular (not packed) */
27276444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
27286444a565SStefano Zampini       PetscInt i, j, N = mumps->id.size_schur;
27296444a565SStefano Zampini       for (i = 0; i < N; i++) {
27306444a565SStefano Zampini         for (j = i; j < N; j++) {
27316444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
27326444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j];
27336444a565SStefano Zampini #else
27343e5b40d0SPierre Jolivet           PetscScalar val = mumps->id.schur[i * N + j].r + PETSC_i * mumps->id.schur[i * N + j].i;
27356444a565SStefano Zampini #endif
27363e5b40d0SPierre Jolivet           array[i * N + j] = array[j * N + i] = val;
27376444a565SStefano Zampini         }
27386444a565SStefano Zampini       }
27396444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
27409566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array, mumps->id.schur, mumps->id.size_schur * mumps->id.size_schur));
27416444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
27426444a565SStefano Zampini       PetscInt i, j, N = mumps->id.size_schur;
27436444a565SStefano Zampini       for (i = 0; i < N; i++) {
27446444a565SStefano Zampini         for (j = 0; j < i + 1; j++) {
27456444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
27466444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i * N + j];
27476444a565SStefano Zampini #else
27483e5b40d0SPierre Jolivet           PetscScalar val = mumps->id.schur[i * N + j].r + PETSC_i * mumps->id.schur[i * N + j].i;
27496444a565SStefano Zampini #endif
27503e5b40d0SPierre Jolivet           array[i * N + j] = array[j * N + i] = val;
27516444a565SStefano Zampini         }
27526444a565SStefano Zampini       }
27536444a565SStefano Zampini     }
27546444a565SStefano Zampini   }
27559566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(St, &array));
27566444a565SStefano Zampini   *S = St;
27573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27586444a565SStefano Zampini }
27596444a565SStefano Zampini 
276066976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F, PetscInt icntl, PetscInt ival)
2761d71ae5a4SJacob Faibussowitsch {
2762e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
27635ccb76cbSHong Zhang 
27645ccb76cbSHong Zhang   PetscFunctionBegin;
2765413bcc21SPierre Jolivet   if (mumps->id.job == JOB_NULL) {                                            /* need to cache icntl and ival since PetscMUMPS_c() has never been called */
27666497c311SBarry Smith     PetscMUMPSInt i, nICNTL_pre = mumps->ICNTL_pre ? mumps->ICNTL_pre[0] : 0; /* number of already cached ICNTL */
27679371c9d4SSatish Balay     for (i = 0; i < nICNTL_pre; ++i)
27689371c9d4SSatish Balay       if (mumps->ICNTL_pre[1 + 2 * i] == icntl) break; /* is this ICNTL already cached? */
2769413bcc21SPierre Jolivet     if (i == nICNTL_pre) {                             /* not already cached */
2770413bcc21SPierre Jolivet       if (i > 0) PetscCall(PetscRealloc(sizeof(PetscMUMPSInt) * (2 * nICNTL_pre + 3), &mumps->ICNTL_pre));
2771413bcc21SPierre Jolivet       else PetscCall(PetscCalloc(sizeof(PetscMUMPSInt) * 3, &mumps->ICNTL_pre));
2772413bcc21SPierre Jolivet       mumps->ICNTL_pre[0]++;
2773413bcc21SPierre Jolivet     }
27746497c311SBarry Smith     mumps->ICNTL_pre[1 + 2 * i] = (PetscMUMPSInt)icntl;
2775413bcc21SPierre Jolivet     PetscCall(PetscMUMPSIntCast(ival, mumps->ICNTL_pre + 2 + 2 * i));
2776413bcc21SPierre Jolivet   } else PetscCall(PetscMUMPSIntCast(ival, &mumps->id.ICNTL(icntl)));
27773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27785ccb76cbSHong Zhang }
27795ccb76cbSHong Zhang 
278066976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F, PetscInt icntl, PetscInt *ival)
2781d71ae5a4SJacob Faibussowitsch {
2782e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2783bc6112feSHong Zhang 
2784bc6112feSHong Zhang   PetscFunctionBegin;
278536df9881Sjeremy theler   if (mumps->id.job == JOB_NULL) {
278636df9881Sjeremy theler     PetscInt i, nICNTL_pre = mumps->ICNTL_pre ? mumps->ICNTL_pre[0] : 0;
278736df9881Sjeremy theler     *ival = 0;
278836df9881Sjeremy theler     for (i = 0; i < nICNTL_pre; ++i) {
278936df9881Sjeremy theler       if (mumps->ICNTL_pre[1 + 2 * i] == icntl) *ival = mumps->ICNTL_pre[2 + 2 * i];
279036df9881Sjeremy theler     }
279136df9881Sjeremy theler   } else *ival = mumps->id.ICNTL(icntl);
27923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2793bc6112feSHong Zhang }
2794bc6112feSHong Zhang 
27955ccb76cbSHong Zhang /*@
27961d27aa22SBarry Smith   MatMumpsSetIcntl - Set MUMPS parameter ICNTL() <https://mumps-solver.org/index.php?page=doc>
27975ccb76cbSHong Zhang 
2798c3339decSBarry Smith   Logically Collective
27995ccb76cbSHong Zhang 
28005ccb76cbSHong Zhang   Input Parameters:
28010b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
28025ccb76cbSHong Zhang . icntl - index of MUMPS parameter array ICNTL()
28035ccb76cbSHong Zhang - ival  - value of MUMPS ICNTL(icntl)
28045ccb76cbSHong Zhang 
28053c7db156SBarry Smith   Options Database Key:
2806147403d9SBarry Smith . -mat_mumps_icntl_<icntl> <ival> - change the option numbered icntl to ival
28075ccb76cbSHong Zhang 
28085ccb76cbSHong Zhang   Level: beginner
28095ccb76cbSHong Zhang 
28101cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
28115ccb76cbSHong Zhang @*/
2812d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsSetIcntl(Mat F, PetscInt icntl, PetscInt ival)
2813d71ae5a4SJacob Faibussowitsch {
28145ccb76cbSHong Zhang   PetscFunctionBegin;
28152989dfd4SHong Zhang   PetscValidType(F, 1);
281628b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
28175ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
28185ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F, ival, 3);
2819c92b4f89SPierre Jolivet   PetscCheck((icntl >= 1 && icntl <= 38) || icntl == 48 || icntl == 58, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported ICNTL value %" PetscInt_FMT, icntl);
2820cac4c232SBarry Smith   PetscTryMethod(F, "MatMumpsSetIcntl_C", (Mat, PetscInt, PetscInt), (F, icntl, ival));
28213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28225ccb76cbSHong Zhang }
28235ccb76cbSHong Zhang 
2824a21f80fcSHong Zhang /*@
28251d27aa22SBarry Smith   MatMumpsGetIcntl - Get MUMPS parameter ICNTL() <https://mumps-solver.org/index.php?page=doc>
2826a21f80fcSHong Zhang 
2827c3339decSBarry Smith   Logically Collective
2828a21f80fcSHong Zhang 
2829a21f80fcSHong Zhang   Input Parameters:
28300b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
2831a21f80fcSHong Zhang - icntl - index of MUMPS parameter array ICNTL()
2832a21f80fcSHong Zhang 
2833a21f80fcSHong Zhang   Output Parameter:
2834a21f80fcSHong Zhang . ival - value of MUMPS ICNTL(icntl)
2835a21f80fcSHong Zhang 
2836a21f80fcSHong Zhang   Level: beginner
2837a21f80fcSHong Zhang 
28381cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
2839a21f80fcSHong Zhang @*/
2840d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetIcntl(Mat F, PetscInt icntl, PetscInt *ival)
2841d71ae5a4SJacob Faibussowitsch {
2842bc6112feSHong Zhang   PetscFunctionBegin;
28432989dfd4SHong Zhang   PetscValidType(F, 1);
284428b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2845bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
28464f572ea9SToby Isaac   PetscAssertPointer(ival, 3);
2847c92b4f89SPierre Jolivet   PetscCheck((icntl >= 1 && icntl <= 38) || icntl == 48 || icntl == 58, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported ICNTL value %" PetscInt_FMT, icntl);
2848cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetIcntl_C", (Mat, PetscInt, PetscInt *), (F, icntl, ival));
28493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2850bc6112feSHong Zhang }
2851bc6112feSHong Zhang 
285266976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F, PetscInt icntl, PetscReal val)
2853d71ae5a4SJacob Faibussowitsch {
2854e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
28558928b65cSHong Zhang 
28568928b65cSHong Zhang   PetscFunctionBegin;
2857413bcc21SPierre Jolivet   if (mumps->id.job == JOB_NULL) {
2858413bcc21SPierre Jolivet     PetscInt i, nCNTL_pre = mumps->CNTL_pre ? mumps->CNTL_pre[0] : 0;
28599371c9d4SSatish Balay     for (i = 0; i < nCNTL_pre; ++i)
28609371c9d4SSatish Balay       if (mumps->CNTL_pre[1 + 2 * i] == icntl) break;
2861413bcc21SPierre Jolivet     if (i == nCNTL_pre) {
2862413bcc21SPierre Jolivet       if (i > 0) PetscCall(PetscRealloc(sizeof(PetscReal) * (2 * nCNTL_pre + 3), &mumps->CNTL_pre));
2863413bcc21SPierre Jolivet       else PetscCall(PetscCalloc(sizeof(PetscReal) * 3, &mumps->CNTL_pre));
2864413bcc21SPierre Jolivet       mumps->CNTL_pre[0]++;
2865413bcc21SPierre Jolivet     }
2866413bcc21SPierre Jolivet     mumps->CNTL_pre[1 + 2 * i] = icntl;
2867413bcc21SPierre Jolivet     mumps->CNTL_pre[2 + 2 * i] = val;
2868413bcc21SPierre Jolivet   } else mumps->id.CNTL(icntl) = val;
28693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28708928b65cSHong Zhang }
28718928b65cSHong Zhang 
287266976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F, PetscInt icntl, PetscReal *val)
2873d71ae5a4SJacob Faibussowitsch {
2874e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2875bc6112feSHong Zhang 
2876bc6112feSHong Zhang   PetscFunctionBegin;
287736df9881Sjeremy theler   if (mumps->id.job == JOB_NULL) {
287836df9881Sjeremy theler     PetscInt i, nCNTL_pre = mumps->CNTL_pre ? mumps->CNTL_pre[0] : 0;
287936df9881Sjeremy theler     *val = 0.0;
288036df9881Sjeremy theler     for (i = 0; i < nCNTL_pre; ++i) {
288136df9881Sjeremy theler       if (mumps->CNTL_pre[1 + 2 * i] == icntl) *val = mumps->CNTL_pre[2 + 2 * i];
288236df9881Sjeremy theler     }
288336df9881Sjeremy theler   } else *val = mumps->id.CNTL(icntl);
28843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2885bc6112feSHong Zhang }
2886bc6112feSHong Zhang 
28878928b65cSHong Zhang /*@
28881d27aa22SBarry Smith   MatMumpsSetCntl - Set MUMPS parameter CNTL() <https://mumps-solver.org/index.php?page=doc>
28898928b65cSHong Zhang 
2890c3339decSBarry Smith   Logically Collective
28918928b65cSHong Zhang 
28928928b65cSHong Zhang   Input Parameters:
28930b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
28948928b65cSHong Zhang . icntl - index of MUMPS parameter array CNTL()
28958928b65cSHong Zhang - val   - value of MUMPS CNTL(icntl)
28968928b65cSHong Zhang 
28973c7db156SBarry Smith   Options Database Key:
2898147403d9SBarry Smith . -mat_mumps_cntl_<icntl> <val> - change the option numbered icntl to ival
28998928b65cSHong Zhang 
29008928b65cSHong Zhang   Level: beginner
29018928b65cSHong Zhang 
29021cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
29038928b65cSHong Zhang @*/
2904d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsSetCntl(Mat F, PetscInt icntl, PetscReal val)
2905d71ae5a4SJacob Faibussowitsch {
29068928b65cSHong Zhang   PetscFunctionBegin;
29072989dfd4SHong Zhang   PetscValidType(F, 1);
290828b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
29098928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
2910bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F, val, 3);
2911413bcc21SPierre Jolivet   PetscCheck(icntl >= 1 && icntl <= 7, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported CNTL value %" PetscInt_FMT, icntl);
2912cac4c232SBarry Smith   PetscTryMethod(F, "MatMumpsSetCntl_C", (Mat, PetscInt, PetscReal), (F, icntl, val));
29133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29148928b65cSHong Zhang }
29158928b65cSHong Zhang 
2916a21f80fcSHong Zhang /*@
29171d27aa22SBarry Smith   MatMumpsGetCntl - Get MUMPS parameter CNTL() <https://mumps-solver.org/index.php?page=doc>
2918a21f80fcSHong Zhang 
2919c3339decSBarry Smith   Logically Collective
2920a21f80fcSHong Zhang 
2921a21f80fcSHong Zhang   Input Parameters:
29220b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
2923a21f80fcSHong Zhang - icntl - index of MUMPS parameter array CNTL()
2924a21f80fcSHong Zhang 
2925a21f80fcSHong Zhang   Output Parameter:
2926a21f80fcSHong Zhang . val - value of MUMPS CNTL(icntl)
2927a21f80fcSHong Zhang 
2928a21f80fcSHong Zhang   Level: beginner
2929a21f80fcSHong Zhang 
29301cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
2931a21f80fcSHong Zhang @*/
2932d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetCntl(Mat F, PetscInt icntl, PetscReal *val)
2933d71ae5a4SJacob Faibussowitsch {
2934bc6112feSHong Zhang   PetscFunctionBegin;
29352989dfd4SHong Zhang   PetscValidType(F, 1);
293628b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
2937bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F, icntl, 2);
29384f572ea9SToby Isaac   PetscAssertPointer(val, 3);
2939413bcc21SPierre Jolivet   PetscCheck(icntl >= 1 && icntl <= 7, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONG, "Unsupported CNTL value %" PetscInt_FMT, icntl);
2940cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetCntl_C", (Mat, PetscInt, PetscReal *), (F, icntl, val));
29413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2942bc6112feSHong Zhang }
2943bc6112feSHong Zhang 
294466976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F, PetscInt icntl, PetscInt *info)
2945d71ae5a4SJacob Faibussowitsch {
2946e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2947bc6112feSHong Zhang 
2948bc6112feSHong Zhang   PetscFunctionBegin;
2949bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
29503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2951bc6112feSHong Zhang }
2952bc6112feSHong Zhang 
295366976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F, PetscInt icntl, PetscInt *infog)
2954d71ae5a4SJacob Faibussowitsch {
2955e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2956bc6112feSHong Zhang 
2957bc6112feSHong Zhang   PetscFunctionBegin;
2958bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
29593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2960bc6112feSHong Zhang }
2961bc6112feSHong Zhang 
296266976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F, PetscInt icntl, PetscReal *rinfo)
2963d71ae5a4SJacob Faibussowitsch {
2964e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2965bc6112feSHong Zhang 
2966bc6112feSHong Zhang   PetscFunctionBegin;
2967bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
29683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2969bc6112feSHong Zhang }
2970bc6112feSHong Zhang 
297166976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F, PetscInt icntl, PetscReal *rinfog)
2972d71ae5a4SJacob Faibussowitsch {
2973e69c285eSBarry Smith   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
2974bc6112feSHong Zhang 
2975bc6112feSHong Zhang   PetscFunctionBegin;
2976bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
29773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2978bc6112feSHong Zhang }
2979bc6112feSHong Zhang 
298066976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetNullPivots_MUMPS(Mat F, PetscInt *size, PetscInt **array)
29815c0bae8cSAshish Patel {
29825c0bae8cSAshish Patel   Mat_MUMPS *mumps = (Mat_MUMPS *)F->data;
29835c0bae8cSAshish Patel 
29845c0bae8cSAshish Patel   PetscFunctionBegin;
29855c0bae8cSAshish 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");
29865c0bae8cSAshish Patel   *size  = 0;
29875c0bae8cSAshish Patel   *array = NULL;
29885c0bae8cSAshish Patel   if (!mumps->myid) {
29895c0bae8cSAshish Patel     *size = mumps->id.INFOG(28);
29905c0bae8cSAshish Patel     PetscCall(PetscMalloc1(*size, array));
29915c0bae8cSAshish Patel     for (int i = 0; i < *size; i++) (*array)[i] = mumps->id.pivnul_list[i] - 1;
29925c0bae8cSAshish Patel   }
29935c0bae8cSAshish Patel   PetscFunctionReturn(PETSC_SUCCESS);
29945c0bae8cSAshish Patel }
29955c0bae8cSAshish Patel 
299666976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F, Mat spRHS)
2997d71ae5a4SJacob Faibussowitsch {
29980e6b8875SHong Zhang   Mat          Bt = NULL, Btseq = NULL;
29990e6b8875SHong Zhang   PetscBool    flg;
3000bb599dfdSHong Zhang   Mat_MUMPS   *mumps = (Mat_MUMPS *)F->data;
3001bb599dfdSHong Zhang   PetscScalar *aa;
3002f410b75aSHong Zhang   PetscInt     spnr, *ia, *ja, M, nrhs;
3003bb599dfdSHong Zhang 
3004bb599dfdSHong Zhang   PetscFunctionBegin;
30054f572ea9SToby Isaac   PetscAssertPointer(spRHS, 2);
3006013e2dc7SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)spRHS, MATTRANSPOSEVIRTUAL, &flg));
30070e6b8875SHong Zhang   if (flg) {
30089566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(spRHS, &Bt));
3009013e2dc7SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)spRHS), PETSC_ERR_ARG_WRONG, "Matrix spRHS must be type MATTRANSPOSEVIRTUAL matrix");
3010bb599dfdSHong Zhang 
30119566063dSJacob Faibussowitsch   PetscCall(MatMumpsSetIcntl(F, 30, 1));
3012bb599dfdSHong Zhang 
30132d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
30140e6b8875SHong Zhang     Mat_MPIAIJ *b = (Mat_MPIAIJ *)Bt->data;
30150e6b8875SHong Zhang     Btseq         = b->A;
30160e6b8875SHong Zhang   } else {
30170e6b8875SHong Zhang     Btseq = Bt;
30180e6b8875SHong Zhang   }
30190e6b8875SHong Zhang 
30209566063dSJacob Faibussowitsch   PetscCall(MatGetSize(spRHS, &M, &nrhs));
30216497c311SBarry Smith   mumps->id.nrhs = (PetscMUMPSInt)nrhs;
30226497c311SBarry Smith   PetscCall(PetscMUMPSIntCast(M, &mumps->id.lrhs));
3023f410b75aSHong Zhang   mumps->id.rhs = NULL;
3024f410b75aSHong Zhang 
3025e3f2db6aSHong Zhang   if (!mumps->myid) {
30269566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArray(Btseq, &aa));
30279566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(Btseq, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
302828b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
30299566063dSJacob Faibussowitsch     PetscCall(PetscMUMPSIntCSRCast(mumps, spnr, ia, ja, &mumps->id.irhs_ptr, &mumps->id.irhs_sparse, &mumps->id.nz_rhs));
3030bb599dfdSHong Zhang     mumps->id.rhs_sparse = (MumpsScalar *)aa;
3031e3f2db6aSHong Zhang   } else {
3032e3f2db6aSHong Zhang     mumps->id.irhs_ptr    = NULL;
3033e3f2db6aSHong Zhang     mumps->id.irhs_sparse = NULL;
3034e3f2db6aSHong Zhang     mumps->id.nz_rhs      = 0;
3035e3f2db6aSHong Zhang     mumps->id.rhs_sparse  = NULL;
3036e3f2db6aSHong Zhang   }
3037bb599dfdSHong Zhang   mumps->id.ICNTL(20) = 1; /* rhs is sparse */
3038e3f2db6aSHong Zhang   mumps->id.ICNTL(21) = 0; /* solution is in assembled centralized format */
3039bb599dfdSHong Zhang 
3040bb599dfdSHong Zhang   /* solve phase */
3041bb599dfdSHong Zhang   mumps->id.job = JOB_SOLVE;
30423ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
30439261f6e4SBarry Smith   PetscCheck(mumps->id.INFOG(1) >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "MUMPS error in solve: INFOG(1)=%d INFO(2)=%d " MUMPS_MANUALS, mumps->id.INFOG(1), mumps->id.INFO(2));
304414267174SHong Zhang 
3045e3f2db6aSHong Zhang   if (!mumps->myid) {
30469566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJRestoreArray(Btseq, &aa));
30479566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(Btseq, 1, PETSC_FALSE, PETSC_FALSE, &spnr, (const PetscInt **)&ia, (const PetscInt **)&ja, &flg));
304828b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot get IJ structure");
3049e3f2db6aSHong Zhang   }
30503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3051bb599dfdSHong Zhang }
3052bb599dfdSHong Zhang 
3053bb599dfdSHong Zhang /*@
30541d27aa22SBarry Smith   MatMumpsGetInverse - Get user-specified set of entries in inverse of `A` <https://mumps-solver.org/index.php?page=doc>
3055bb599dfdSHong Zhang 
3056c3339decSBarry Smith   Logically Collective
3057bb599dfdSHong Zhang 
305820f4b53cSBarry Smith   Input Parameter:
30590b4b7b1cSBarry Smith . F - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
3060bb599dfdSHong Zhang 
3061bb599dfdSHong Zhang   Output Parameter:
306220f4b53cSBarry Smith . spRHS - sequential sparse matrix in `MATTRANSPOSEVIRTUAL` format with requested entries of inverse of `A`
3063bb599dfdSHong Zhang 
3064bb599dfdSHong Zhang   Level: beginner
3065bb599dfdSHong Zhang 
30661cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatCreateTranspose()`
3067bb599dfdSHong Zhang @*/
3068d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInverse(Mat F, Mat spRHS)
3069d71ae5a4SJacob Faibussowitsch {
3070bb599dfdSHong Zhang   PetscFunctionBegin;
3071bb599dfdSHong Zhang   PetscValidType(F, 1);
307228b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
3073cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInverse_C", (Mat, Mat), (F, spRHS));
30743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3075bb599dfdSHong Zhang }
3076bb599dfdSHong Zhang 
307766976f2fSJacob Faibussowitsch static PetscErrorCode MatMumpsGetInverseTranspose_MUMPS(Mat F, Mat spRHST)
3078d71ae5a4SJacob Faibussowitsch {
30790e6b8875SHong Zhang   Mat spRHS;
30800e6b8875SHong Zhang 
30810e6b8875SHong Zhang   PetscFunctionBegin;
30829566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(spRHST, &spRHS));
30839566063dSJacob Faibussowitsch   PetscCall(MatMumpsGetInverse_MUMPS(F, spRHS));
30849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&spRHS));
30853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30860e6b8875SHong Zhang }
30870e6b8875SHong Zhang 
30880e6b8875SHong Zhang /*@
30891d27aa22SBarry Smith   MatMumpsGetInverseTranspose - Get user-specified set of entries in inverse of matrix $A^T $ <https://mumps-solver.org/index.php?page=doc>
30900e6b8875SHong Zhang 
3091c3339decSBarry Smith   Logically Collective
30920e6b8875SHong Zhang 
309320f4b53cSBarry Smith   Input Parameter:
30940b4b7b1cSBarry Smith . F - the factored matrix of A obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
30950e6b8875SHong Zhang 
30960e6b8875SHong Zhang   Output Parameter:
309720f4b53cSBarry Smith . spRHST - sequential sparse matrix in `MATAIJ` format containing the requested entries of inverse of `A`^T
30980e6b8875SHong Zhang 
30990e6b8875SHong Zhang   Level: beginner
31000e6b8875SHong Zhang 
31011cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatCreateTranspose()`, `MatMumpsGetInverse()`
31020e6b8875SHong Zhang @*/
3103d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInverseTranspose(Mat F, Mat spRHST)
3104d71ae5a4SJacob Faibussowitsch {
31050e6b8875SHong Zhang   PetscBool flg;
31060e6b8875SHong Zhang 
31070e6b8875SHong Zhang   PetscFunctionBegin;
31080e6b8875SHong Zhang   PetscValidType(F, 1);
310928b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
31109566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)spRHST, &flg, MATSEQAIJ, MATMPIAIJ, NULL));
311128b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)spRHST), PETSC_ERR_ARG_WRONG, "Matrix spRHST must be MATAIJ matrix");
31120e6b8875SHong Zhang 
3113cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInverseTranspose_C", (Mat, Mat), (F, spRHST));
31143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31150e6b8875SHong Zhang }
31160e6b8875SHong Zhang 
3117a21f80fcSHong Zhang /*@
31181d27aa22SBarry Smith   MatMumpsGetInfo - Get MUMPS parameter INFO() <https://mumps-solver.org/index.php?page=doc>
3119a21f80fcSHong Zhang 
3120c3339decSBarry Smith   Logically Collective
3121a21f80fcSHong Zhang 
3122a21f80fcSHong Zhang   Input Parameters:
31230b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
3124a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFO()
3125a21f80fcSHong Zhang 
3126a21f80fcSHong Zhang   Output Parameter:
3127a21f80fcSHong Zhang . ival - value of MUMPS INFO(icntl)
3128a21f80fcSHong Zhang 
3129a21f80fcSHong Zhang   Level: beginner
3130a21f80fcSHong Zhang 
31311cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
3132a21f80fcSHong Zhang @*/
3133d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInfo(Mat F, PetscInt icntl, PetscInt *ival)
3134d71ae5a4SJacob Faibussowitsch {
3135bc6112feSHong Zhang   PetscFunctionBegin;
31362989dfd4SHong Zhang   PetscValidType(F, 1);
313728b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
31384f572ea9SToby Isaac   PetscAssertPointer(ival, 3);
3139cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInfo_C", (Mat, PetscInt, PetscInt *), (F, icntl, ival));
31403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3141bc6112feSHong Zhang }
3142bc6112feSHong Zhang 
3143a21f80fcSHong Zhang /*@
31441d27aa22SBarry Smith   MatMumpsGetInfog - Get MUMPS parameter INFOG() <https://mumps-solver.org/index.php?page=doc>
3145a21f80fcSHong Zhang 
3146c3339decSBarry Smith   Logically Collective
3147a21f80fcSHong Zhang 
3148a21f80fcSHong Zhang   Input Parameters:
31490b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
3150a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFOG()
3151a21f80fcSHong Zhang 
3152a21f80fcSHong Zhang   Output Parameter:
3153a21f80fcSHong Zhang . ival - value of MUMPS INFOG(icntl)
3154a21f80fcSHong Zhang 
3155a21f80fcSHong Zhang   Level: beginner
3156a21f80fcSHong Zhang 
31571cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`
3158a21f80fcSHong Zhang @*/
3159d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetInfog(Mat F, PetscInt icntl, PetscInt *ival)
3160d71ae5a4SJacob Faibussowitsch {
3161bc6112feSHong Zhang   PetscFunctionBegin;
31622989dfd4SHong Zhang   PetscValidType(F, 1);
316328b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
31644f572ea9SToby Isaac   PetscAssertPointer(ival, 3);
3165cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetInfog_C", (Mat, PetscInt, PetscInt *), (F, icntl, ival));
31663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3167bc6112feSHong Zhang }
3168bc6112feSHong Zhang 
3169a21f80fcSHong Zhang /*@
31701d27aa22SBarry Smith   MatMumpsGetRinfo - Get MUMPS parameter RINFO() <https://mumps-solver.org/index.php?page=doc>
3171a21f80fcSHong Zhang 
3172c3339decSBarry Smith   Logically Collective
3173a21f80fcSHong Zhang 
3174a21f80fcSHong Zhang   Input Parameters:
31750b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
3176a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFO()
3177a21f80fcSHong Zhang 
3178a21f80fcSHong Zhang   Output Parameter:
3179a21f80fcSHong Zhang . val - value of MUMPS RINFO(icntl)
3180a21f80fcSHong Zhang 
3181a21f80fcSHong Zhang   Level: beginner
3182a21f80fcSHong Zhang 
31831cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfog()`
3184a21f80fcSHong Zhang @*/
3185d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetRinfo(Mat F, PetscInt icntl, PetscReal *val)
3186d71ae5a4SJacob Faibussowitsch {
3187bc6112feSHong Zhang   PetscFunctionBegin;
31882989dfd4SHong Zhang   PetscValidType(F, 1);
318928b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
31904f572ea9SToby Isaac   PetscAssertPointer(val, 3);
3191cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetRinfo_C", (Mat, PetscInt, PetscReal *), (F, icntl, val));
31923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3193bc6112feSHong Zhang }
3194bc6112feSHong Zhang 
3195a21f80fcSHong Zhang /*@
31961d27aa22SBarry Smith   MatMumpsGetRinfog - Get MUMPS parameter RINFOG() <https://mumps-solver.org/index.php?page=doc>
3197a21f80fcSHong Zhang 
3198c3339decSBarry Smith   Logically Collective
3199a21f80fcSHong Zhang 
3200a21f80fcSHong Zhang   Input Parameters:
32010b4b7b1cSBarry Smith + F     - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
3202a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFOG()
3203a21f80fcSHong Zhang 
3204a21f80fcSHong Zhang   Output Parameter:
3205a21f80fcSHong Zhang . val - value of MUMPS RINFOG(icntl)
3206a21f80fcSHong Zhang 
3207a21f80fcSHong Zhang   Level: beginner
3208a21f80fcSHong Zhang 
32091cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`
3210a21f80fcSHong Zhang @*/
3211d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMumpsGetRinfog(Mat F, PetscInt icntl, PetscReal *val)
3212d71ae5a4SJacob Faibussowitsch {
3213bc6112feSHong Zhang   PetscFunctionBegin;
32142989dfd4SHong Zhang   PetscValidType(F, 1);
321528b400f6SJacob Faibussowitsch   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
32164f572ea9SToby Isaac   PetscAssertPointer(val, 3);
3217cac4c232SBarry Smith   PetscUseMethod(F, "MatMumpsGetRinfog_C", (Mat, PetscInt, PetscReal *), (F, icntl, val));
32183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3219bc6112feSHong Zhang }
3220bc6112feSHong Zhang 
32215c0bae8cSAshish Patel /*@
32221d27aa22SBarry Smith   MatMumpsGetNullPivots - Get MUMPS parameter PIVNUL_LIST() <https://mumps-solver.org/index.php?page=doc>
32235c0bae8cSAshish Patel 
32245c0bae8cSAshish Patel   Logically Collective
32255c0bae8cSAshish Patel 
32265c0bae8cSAshish Patel   Input Parameter:
32270b4b7b1cSBarry Smith . F - the factored matrix obtained by calling `MatGetFactor()` with a `MatSolverType` of `MATSOLVERMUMPS` and a `MatFactorType` of `MAT_FACTOR_LU` or `MAT_FACTOR_CHOLESKY`
32285c0bae8cSAshish Patel 
32295c0bae8cSAshish Patel   Output Parameters:
32300b4b7b1cSBarry Smith + size  - local size of the array. The size of the array is non-zero only on MPI rank 0
32315c0bae8cSAshish 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
32325c0bae8cSAshish Patel           for freeing this array.
32335c0bae8cSAshish Patel 
32345c0bae8cSAshish Patel   Level: beginner
32355c0bae8cSAshish Patel 
32361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`
32375c0bae8cSAshish Patel @*/
32385c0bae8cSAshish Patel PetscErrorCode MatMumpsGetNullPivots(Mat F, PetscInt *size, PetscInt **array)
32395c0bae8cSAshish Patel {
32405c0bae8cSAshish Patel   PetscFunctionBegin;
32415c0bae8cSAshish Patel   PetscValidType(F, 1);
32425c0bae8cSAshish Patel   PetscCheck(F->factortype, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
32434f572ea9SToby Isaac   PetscAssertPointer(size, 2);
32444f572ea9SToby Isaac   PetscAssertPointer(array, 3);
32455c0bae8cSAshish Patel   PetscUseMethod(F, "MatMumpsGetNullPivots_C", (Mat, PetscInt *, PetscInt **), (F, size, array));
32465c0bae8cSAshish Patel   PetscFunctionReturn(PETSC_SUCCESS);
32475c0bae8cSAshish Patel }
32485c0bae8cSAshish Patel 
324924b6179bSKris Buschelman /*MC
32502692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
32510b4b7b1cSBarry Smith   MPI distributed and sequential matrices via the external package MUMPS <https://mumps-solver.org/index.php?page=doc>
325224b6179bSKris Buschelman 
325311a5261eSBarry Smith   Works with `MATAIJ` and `MATSBAIJ` matrices
325424b6179bSKris Buschelman 
3255c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS
3256c2b89b5dSBarry Smith 
32572ef1f0ffSBarry 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.
32582ef1f0ffSBarry Smith   See details below.
3259217d3b1eSJunchao Zhang 
32602ef1f0ffSBarry Smith   Use `-pc_type cholesky` or `lu` `-pc_factor_mat_solver_type mumps` to use this direct solver
3261c2b89b5dSBarry Smith 
326224b6179bSKris Buschelman   Options Database Keys:
32634422a9fcSPatrick Sanan +  -mat_mumps_icntl_1  - ICNTL(1): output stream for error messages
32644422a9fcSPatrick Sanan .  -mat_mumps_icntl_2  - ICNTL(2): output stream for diagnostic printing, statistics, and warning
32654422a9fcSPatrick Sanan .  -mat_mumps_icntl_3  - ICNTL(3): output stream for global information, collected on the host
32664422a9fcSPatrick Sanan .  -mat_mumps_icntl_4  - ICNTL(4): level of printing (0 to 4)
32674422a9fcSPatrick Sanan .  -mat_mumps_icntl_6  - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
3268b53c1a7fSBarry 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
3269b53c1a7fSBarry Smith                           Use -pc_factor_mat_ordering_type <type> to have PETSc perform the ordering (sequential only)
32704422a9fcSPatrick Sanan .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
32714422a9fcSPatrick Sanan .  -mat_mumps_icntl_10 - ICNTL(10): max num of refinements
32724422a9fcSPatrick Sanan .  -mat_mumps_icntl_11 - ICNTL(11): statistics related to an error analysis (via -ksp_view)
32734422a9fcSPatrick Sanan .  -mat_mumps_icntl_12 - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
32744422a9fcSPatrick Sanan .  -mat_mumps_icntl_13 - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
32754422a9fcSPatrick Sanan .  -mat_mumps_icntl_14 - ICNTL(14): percentage increase in the estimated working space
327645e3843bSPierre Jolivet .  -mat_mumps_icntl_15 - ICNTL(15): compression of the input matrix resulting from a block format
32774422a9fcSPatrick Sanan .  -mat_mumps_icntl_19 - ICNTL(19): computes the Schur complement
327825aac85cSJunchao Zhang .  -mat_mumps_icntl_20 - ICNTL(20): give MUMPS centralized (0) or distributed (10) dense RHS
32794422a9fcSPatrick Sanan .  -mat_mumps_icntl_22 - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
32804422a9fcSPatrick Sanan .  -mat_mumps_icntl_23 - ICNTL(23): max size of the working memory (MB) that can allocate per processor
32814422a9fcSPatrick Sanan .  -mat_mumps_icntl_24 - ICNTL(24): detection of null pivot rows (0 or 1)
32824422a9fcSPatrick Sanan .  -mat_mumps_icntl_25 - ICNTL(25): compute a solution of a deficient matrix and a null space basis
32834422a9fcSPatrick Sanan .  -mat_mumps_icntl_26 - ICNTL(26): drives the solution phase if a Schur complement matrix
3284fa6fd9d0SPierre Jolivet .  -mat_mumps_icntl_28 - ICNTL(28): use 1 for sequential analysis and ICNTL(7) ordering, or 2 for parallel analysis and ICNTL(29) ordering
32854422a9fcSPatrick Sanan .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
32864422a9fcSPatrick Sanan .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
32874422a9fcSPatrick Sanan .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
32884422a9fcSPatrick Sanan .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
3289a0e18203SThibaut Appel .  -mat_mumps_icntl_35 - ICNTL(35): level of activation of BLR (Block Low-Rank) feature
3290a0e18203SThibaut Appel .  -mat_mumps_icntl_36 - ICNTL(36): controls the choice of BLR factorization variant
3291*50ea2040Saszaboa .  -mat_mumps_icntl_37 - ICNTL(37): compression of the contribution blocks (CB)
3292a0e18203SThibaut Appel .  -mat_mumps_icntl_38 - ICNTL(38): sets the estimated compression rate of LU factors with BLR
3293c92b4f89SPierre Jolivet .  -mat_mumps_icntl_48 - ICNTL(48): multithreading with tree parallelism
3294146931dbSPierre Jolivet .  -mat_mumps_icntl_58 - ICNTL(58): options for symbolic factorization
32954422a9fcSPatrick Sanan .  -mat_mumps_cntl_1   - CNTL(1): relative pivoting threshold
32964422a9fcSPatrick Sanan .  -mat_mumps_cntl_2   - CNTL(2): stopping criterion of refinement
32974422a9fcSPatrick Sanan .  -mat_mumps_cntl_3   - CNTL(3): absolute pivoting threshold
32984422a9fcSPatrick Sanan .  -mat_mumps_cntl_4   - CNTL(4): value for static pivoting
3299217d3b1eSJunchao Zhang .  -mat_mumps_cntl_5   - CNTL(5): fixation for null pivots
3300a0e18203SThibaut Appel .  -mat_mumps_cntl_7   - CNTL(7): precision of the dropping parameter used during BLR factorization
3301217d3b1eSJunchao 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.
3302217d3b1eSJunchao Zhang                                     Default might be the number of cores per CPU package (socket) as reported by hwloc and suggested by the MUMPS manual.
330324b6179bSKris Buschelman 
330424b6179bSKris Buschelman   Level: beginner
330524b6179bSKris Buschelman 
330695452b02SPatrick Sanan   Notes:
33071d27aa22SBarry 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
33082ef1f0ffSBarry Smith   error if the matrix is Hermitian.
330938548759SBarry Smith 
331026cc229bSBarry 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
331126cc229bSBarry Smith   `MatSetOptionsPrefixFactor()` on the matrix from which the factor was obtained or `MatSetOptionsPrefix()` on the factor matrix.
331226cc229bSBarry Smith 
33132ef1f0ffSBarry 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
33142ef1f0ffSBarry Smith   the failure with
33152ef1f0ffSBarry Smith .vb
33162ef1f0ffSBarry Smith           KSPGetPC(ksp,&pc);
33172ef1f0ffSBarry Smith           PCFactorGetMatrix(pc,&mat);
33182ef1f0ffSBarry Smith           MatMumpsGetInfo(mat,....);
33192ef1f0ffSBarry Smith           MatMumpsGetInfog(mat,....); etc.
33202ef1f0ffSBarry Smith .ve
33212ef1f0ffSBarry Smith   Or run with `-ksp_error_if_not_converged` and the program will be stopped and the information printed in the error message.
33229fc87aa7SBarry Smith 
3323a5399872SJunchao Zhang   MUMPS provides 64-bit integer support in two build modes:
3324a5399872SJunchao Zhang   full 64-bit: here MUMPS is built with C preprocessing flag -DINTSIZE64 and Fortran compiler option -i8, -fdefault-integer-8 or equivalent, and
3325a5399872SJunchao 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).
33268fcaa860SBarry Smith 
3327a5399872SJunchao 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,
3328a5399872SJunchao 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
3329a5399872SJunchao 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
3330a5399872SJunchao Zhang   integer) build of all dependent libraries MPI, ScaLAPACK, LAPACK and BLAS.
3331a5399872SJunchao Zhang 
3332a5399872SJunchao 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.
3333a5399872SJunchao Zhang 
3334a5399872SJunchao Zhang   Two modes to run MUMPS/PETSc with OpenMP
33352ef1f0ffSBarry Smith .vb
33360b4b7b1cSBarry Smith    Set `OMP_NUM_THREADS` and run with fewer MPI ranks than cores. For example, if you want to have 16 OpenMP
33370b4b7b1cSBarry Smith    threads per rank, then you may use "export `OMP_NUM_THREADS` = 16 && mpirun -n 4 ./test".
33382ef1f0ffSBarry Smith .ve
33398fcaa860SBarry Smith 
33402ef1f0ffSBarry Smith .vb
33410b4b7b1cSBarry Smith    `-mat_mumps_use_omp_threads` [m] and run your code with as many MPI ranks as the number of cores. For example,
33422ef1f0ffSBarry 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"
33432ef1f0ffSBarry Smith .ve
33448fcaa860SBarry Smith 
33458fcaa860SBarry Smith    To run MUMPS in MPI+OpenMP hybrid mode (i.e., enable multithreading in MUMPS), but still run the non-MUMPS part
33462ef1f0ffSBarry 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`
33472ef1f0ffSBarry 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
33488fcaa860SBarry Smith    libraries, to really get performance, you should have multithreaded BLAS libraries such as Intel MKL, AMD ACML, Cray libSci or OpenBLAS
33490b4b7b1cSBarry Smith    (PETSc will automatically try to utilized a threaded BLAS if `--with-openmp` is provided).
3350217d3b1eSJunchao Zhang 
33518fcaa860SBarry 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
3352217d3b1eSJunchao Zhang    processes on each compute node. Listing the processes in rank ascending order, we split processes on a node into consecutive groups of
3353217d3b1eSJunchao 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
3354217d3b1eSJunchao 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
3355217d3b1eSJunchao 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.
3356217d3b1eSJunchao 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,
3357217d3b1eSJunchao 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
33580b4b7b1cSBarry Smith    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
3359217d3b1eSJunchao 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
33600b4b7b1cSBarry Smith    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.
33610b4b7b1cSBarry 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
3362217d3b1eSJunchao Zhang    examine the mapping result.
3363217d3b1eSJunchao Zhang 
336411a5261eSBarry 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,
336511a5261eSBarry 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
336611a5261eSBarry Smith    calls `omp_set_num_threads`(m) internally before calling MUMPS.
3367217d3b1eSJunchao Zhang 
33681d27aa22SBarry Smith    See {cite}`heroux2011bi` and {cite}`gutierrez2017accommodating`
3369217d3b1eSJunchao Zhang 
33701cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatMumpsSetIcntl()`, `MatMumpsGetIcntl()`, `MatMumpsSetCntl()`, `MatMumpsGetCntl()`, `MatMumpsGetInfo()`, `MatMumpsGetInfog()`, `MatMumpsGetRinfo()`, `MatMumpsGetRinfog()`, `KSPGetPC()`, `PCFactorGetMatrix()`
337124b6179bSKris Buschelman M*/
337224b6179bSKris Buschelman 
3373d2a308c1SPierre Jolivet static PetscErrorCode MatFactorGetSolverType_mumps(PETSC_UNUSED Mat A, MatSolverType *type)
3374d71ae5a4SJacob Faibussowitsch {
337535bd34faSBarry Smith   PetscFunctionBegin;
33762692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
33773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337835bd34faSBarry Smith }
337935bd34faSBarry Smith 
3380bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
3381d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_aij_mumps(Mat A, MatFactorType ftype, Mat *F)
3382d71ae5a4SJacob Faibussowitsch {
33832877fffaSHong Zhang   Mat         B;
33842877fffaSHong Zhang   Mat_MUMPS  *mumps;
33854b9405b2SPierre Jolivet   PetscBool   isSeqAIJ, isDiag, isDense;
33862c7c0729SBarry Smith   PetscMPIInt size;
33872877fffaSHong Zhang 
33882877fffaSHong Zhang   PetscFunctionBegin;
3389eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
339003e5aca4SStefano Zampini   if (ftype == MAT_FACTOR_CHOLESKY && A->hermitian == PETSC_BOOL3_TRUE && A->symmetric != PETSC_BOOL3_TRUE) {
339103e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Hermitian MAT_FACTOR_CHOLESKY is not supported. Use MAT_FACTOR_LU instead.\n"));
339203e5aca4SStefano Zampini     *F = NULL;
339303e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
339403e5aca4SStefano Zampini   }
3395eb1ec7c1SStefano Zampini #endif
33962877fffaSHong Zhang   /* Create the factorization matrix */
33979566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeqAIJ));
3398c3e1b152SPierre Jolivet   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATDIAGONAL, &isDiag));
33994b9405b2SPierre Jolivet   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &isDense, MATSEQDENSE, MATMPIDENSE, NULL));
34009566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
34019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
3402d2a308c1SPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &((PetscObject)B)->type_name));
34039566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
34042877fffaSHong Zhang 
34054dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
34062205254eSKarl Rupp 
34072877fffaSHong Zhang   B->ops->view    = MatView_MUMPS;
340835bd34faSBarry Smith   B->ops->getinfo = MatGetInfo_MUMPS;
34092205254eSKarl Rupp 
34109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
34119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
34129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
34139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
34149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
34159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
34169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
34179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
34189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
34199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
34209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
34215c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
34229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
34239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
34246444a565SStefano Zampini 
3425450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
3426450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
3427d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
3428bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
3429c3e1b152SPierre Jolivet     else if (isDiag) mumps->ConvertToTriples = MatConvertToTriples_diagonal_xaij;
34304b9405b2SPierre Jolivet     else if (isDense) mumps->ConvertToTriples = MatConvertToTriples_dense_xaij;
3431bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
34329566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
3433746480a1SHong Zhang     mumps->sym = 0;
3434dcd589f8SShri Abhyankar   } else {
343567877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
3436450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
3437bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
3438c3e1b152SPierre Jolivet     else if (isDiag) mumps->ConvertToTriples = MatConvertToTriples_diagonal_xaij;
34394b9405b2SPierre Jolivet     else if (isDense) mumps->ConvertToTriples = MatConvertToTriples_dense_xaij;
3440bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
34419566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
344259ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
344359ac8732SStefano Zampini     mumps->sym = 2;
344459ac8732SStefano Zampini #else
3445b94d7dedSBarry Smith     if (A->spd == PETSC_BOOL3_TRUE) mumps->sym = 1;
34466fdc2a6dSBarry Smith     else mumps->sym = 2;
344759ac8732SStefano Zampini #endif
3448450b117fSShri Abhyankar   }
34492877fffaSHong Zhang 
345000c67f3bSHong Zhang   /* set solvertype */
34519566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
34529566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
34539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
34542c7c0729SBarry Smith   if (size == 1) {
34554ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3456f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
34572c7c0729SBarry Smith   }
34582877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
3459e69c285eSBarry Smith   B->data         = (void *)mumps;
34602205254eSKarl Rupp 
34612877fffaSHong Zhang   *F               = B;
3462413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3463413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3464413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3465d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
34663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34672877fffaSHong Zhang }
34682877fffaSHong Zhang 
3469bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
3470d2a308c1SPierre Jolivet static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A, PETSC_UNUSED MatFactorType ftype, Mat *F)
3471d71ae5a4SJacob Faibussowitsch {
34722877fffaSHong Zhang   Mat         B;
34732877fffaSHong Zhang   Mat_MUMPS  *mumps;
3474ace3abfcSBarry Smith   PetscBool   isSeqSBAIJ;
34752c7c0729SBarry Smith   PetscMPIInt size;
34762877fffaSHong Zhang 
34772877fffaSHong Zhang   PetscFunctionBegin;
3478eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
347903e5aca4SStefano Zampini   if (ftype == MAT_FACTOR_CHOLESKY && A->hermitian == PETSC_BOOL3_TRUE && A->symmetric != PETSC_BOOL3_TRUE) {
348003e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Hermitian MAT_FACTOR_CHOLESKY is not supported. Use MAT_FACTOR_LU instead.\n"));
348103e5aca4SStefano Zampini     *F = NULL;
348203e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
348303e5aca4SStefano Zampini   }
3484eb1ec7c1SStefano Zampini #endif
34859566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
34869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
3487d2a308c1SPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &((PetscObject)B)->type_name));
34889566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3489e69c285eSBarry Smith 
34904dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
34919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSeqSBAIJ));
3492bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
349316ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
3494dcd589f8SShri Abhyankar   } else {
3495bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
3496bccb9932SShri Abhyankar   }
3497bccb9932SShri Abhyankar 
349867877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
3499bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
3500722b6324SPierre Jolivet   B->ops->getinfo                = MatGetInfo_MUMPS;
35012205254eSKarl Rupp 
35029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
35039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
35049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
35059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
35069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
35079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
35089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
35099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
35109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
35119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
35129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
35135c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
35149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
35159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
35162205254eSKarl Rupp 
3517f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
351859ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
351959ac8732SStefano Zampini   mumps->sym = 2;
352059ac8732SStefano Zampini #else
3521b94d7dedSBarry Smith   if (A->spd == PETSC_BOOL3_TRUE) mumps->sym = 1;
35226fdc2a6dSBarry Smith   else mumps->sym = 2;
352359ac8732SStefano Zampini #endif
3524a214ac2aSShri Abhyankar 
352500c67f3bSHong Zhang   /* set solvertype */
35269566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
35279566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
35289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
35292c7c0729SBarry Smith   if (size == 1) {
35304ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3531f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
35322c7c0729SBarry Smith   }
35339566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
3534f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
3535e69c285eSBarry Smith   B->data         = (void *)mumps;
35362205254eSKarl Rupp 
35372877fffaSHong Zhang   *F               = B;
3538413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3539413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3540413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3541d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
35423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35432877fffaSHong Zhang }
354497969023SHong Zhang 
3545d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_baij_mumps(Mat A, MatFactorType ftype, Mat *F)
3546d71ae5a4SJacob Faibussowitsch {
354767877ebaSShri Abhyankar   Mat         B;
354867877ebaSShri Abhyankar   Mat_MUMPS  *mumps;
3549ace3abfcSBarry Smith   PetscBool   isSeqBAIJ;
35502c7c0729SBarry Smith   PetscMPIInt size;
355167877ebaSShri Abhyankar 
355267877ebaSShri Abhyankar   PetscFunctionBegin;
355367877ebaSShri Abhyankar   /* Create the factorization matrix */
35549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &isSeqBAIJ));
35559566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
35569566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
3557d2a308c1SPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &((PetscObject)B)->type_name));
35589566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3559450b117fSShri Abhyankar 
35604dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
3561450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
3562450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
3563450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
3564bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
3565bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
3566746480a1SHong Zhang     mumps->sym = 0;
35679566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
3568546078acSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead");
3569bccb9932SShri Abhyankar 
3570450b117fSShri Abhyankar   B->ops->view    = MatView_MUMPS;
3571722b6324SPierre Jolivet   B->ops->getinfo = MatGetInfo_MUMPS;
35722205254eSKarl Rupp 
35739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
35749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
35759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
35769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
35779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
35789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
35799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
35809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
35819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
35829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
35839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
35845c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
35859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
35869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
3587450b117fSShri Abhyankar 
358800c67f3bSHong Zhang   /* set solvertype */
35899566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
35909566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
35919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
35922c7c0729SBarry Smith   if (size == 1) {
35934ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3594f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
35952c7c0729SBarry Smith   }
35967ee00b23SStefano Zampini   B->ops->destroy = MatDestroy_MUMPS;
35977ee00b23SStefano Zampini   B->data         = (void *)mumps;
35987ee00b23SStefano Zampini 
35997ee00b23SStefano Zampini   *F               = B;
3600413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3601413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3602413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3603d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
36043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36057ee00b23SStefano Zampini }
36067ee00b23SStefano Zampini 
36077ee00b23SStefano Zampini /* MatGetFactor for Seq and MPI SELL matrices */
3608d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_sell_mumps(Mat A, MatFactorType ftype, Mat *F)
3609d71ae5a4SJacob Faibussowitsch {
36107ee00b23SStefano Zampini   Mat         B;
36117ee00b23SStefano Zampini   Mat_MUMPS  *mumps;
36127ee00b23SStefano Zampini   PetscBool   isSeqSELL;
36132c7c0729SBarry Smith   PetscMPIInt size;
36147ee00b23SStefano Zampini 
36157ee00b23SStefano Zampini   PetscFunctionBegin;
36167ee00b23SStefano Zampini   /* Create the factorization matrix */
36179566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSELL, &isSeqSELL));
36189566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
36199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
3620d2a308c1SPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &((PetscObject)B)->type_name));
36219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
36227ee00b23SStefano Zampini 
36234dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mumps));
36247ee00b23SStefano Zampini 
36257ee00b23SStefano Zampini   B->ops->view    = MatView_MUMPS;
36267ee00b23SStefano Zampini   B->ops->getinfo = MatGetInfo_MUMPS;
36277ee00b23SStefano Zampini 
36289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
36299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
36309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
36319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
36329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
36339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
36349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
36359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
36369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
36379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
36389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
36395c0bae8cSAshish Patel   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
36407ee00b23SStefano Zampini 
36417ee00b23SStefano Zampini   if (ftype == MAT_FACTOR_LU) {
36427ee00b23SStefano Zampini     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
36437ee00b23SStefano Zampini     B->factortype            = MAT_FACTOR_LU;
36447ee00b23SStefano Zampini     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
36457ee00b23SStefano Zampini     else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "To be implemented");
36467ee00b23SStefano Zampini     mumps->sym = 0;
36479566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
36487ee00b23SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "To be implemented");
36497ee00b23SStefano Zampini 
36507ee00b23SStefano Zampini   /* set solvertype */
36519566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
36529566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
36539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
36542c7c0729SBarry Smith   if (size == 1) {
36554ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization  */
3656f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
36572c7c0729SBarry Smith   }
3658450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
3659e69c285eSBarry Smith   B->data         = (void *)mumps;
36602205254eSKarl Rupp 
3661450b117fSShri Abhyankar   *F               = B;
3662413bcc21SPierre Jolivet   mumps->id.job    = JOB_NULL;
3663413bcc21SPierre Jolivet   mumps->ICNTL_pre = NULL;
3664413bcc21SPierre Jolivet   mumps->CNTL_pre  = NULL;
3665d47f36abSHong Zhang   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
36663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3667450b117fSShri Abhyankar }
366842c9c57cSBarry Smith 
36699d0448ceSStefano Zampini /* MatGetFactor for MATNEST matrices */
36709d0448ceSStefano Zampini static PetscErrorCode MatGetFactor_nest_mumps(Mat A, MatFactorType ftype, Mat *F)
36719d0448ceSStefano Zampini {
36729d0448ceSStefano Zampini   Mat         B, **mats;
36739d0448ceSStefano Zampini   Mat_MUMPS  *mumps;
36749d0448ceSStefano Zampini   PetscInt    nr, nc;
36759d0448ceSStefano Zampini   PetscMPIInt size;
367603e5aca4SStefano Zampini   PetscBool   flg = PETSC_TRUE;
36779d0448ceSStefano Zampini 
36789d0448ceSStefano Zampini   PetscFunctionBegin;
36799d0448ceSStefano Zampini #if defined(PETSC_USE_COMPLEX)
368003e5aca4SStefano Zampini   if (ftype == MAT_FACTOR_CHOLESKY && A->hermitian == PETSC_BOOL3_TRUE && A->symmetric != PETSC_BOOL3_TRUE) {
368103e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Hermitian MAT_FACTOR_CHOLESKY is not supported. Use MAT_FACTOR_LU instead.\n"));
368203e5aca4SStefano Zampini     *F = NULL;
368303e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
368403e5aca4SStefano Zampini   }
36859d0448ceSStefano Zampini #endif
36869d0448ceSStefano Zampini 
368703e5aca4SStefano Zampini   /* Return if some condition is not satisfied */
368803e5aca4SStefano Zampini   *F = NULL;
36899d0448ceSStefano Zampini   PetscCall(MatNestGetSubMats(A, &nr, &nc, &mats));
36909d0448ceSStefano Zampini   if (ftype == MAT_FACTOR_CHOLESKY) {
36919d0448ceSStefano Zampini     IS       *rows, *cols;
36929d0448ceSStefano Zampini     PetscInt *m, *M;
36939d0448ceSStefano Zampini 
36949d0448ceSStefano Zampini     PetscCheck(nr == nc, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MAT_FACTOR_CHOLESKY not supported for nest sizes %" PetscInt_FMT " != %" PetscInt_FMT ". Use MAT_FACTOR_LU.", nr, nc);
36959d0448ceSStefano Zampini     PetscCall(PetscMalloc2(nr, &rows, nc, &cols));
36969d0448ceSStefano Zampini     PetscCall(MatNestGetISs(A, rows, cols));
36979d0448ceSStefano Zampini     for (PetscInt r = 0; flg && r < nr; r++) PetscCall(ISEqualUnsorted(rows[r], cols[r], &flg));
369803e5aca4SStefano Zampini     if (!flg) {
369903e5aca4SStefano Zampini       PetscCall(PetscFree2(rows, cols));
370003e5aca4SStefano Zampini       PetscCall(PetscInfo(A, "MAT_FACTOR_CHOLESKY not supported for unequal row and column maps. Use MAT_FACTOR_LU.\n"));
370103e5aca4SStefano Zampini       PetscFunctionReturn(PETSC_SUCCESS);
370203e5aca4SStefano Zampini     }
37039d0448ceSStefano Zampini     PetscCall(PetscMalloc2(nr, &m, nr, &M));
37049d0448ceSStefano Zampini     for (PetscInt r = 0; r < nr; r++) PetscCall(ISGetMinMax(rows[r], &m[r], &M[r]));
37059d0448ceSStefano Zampini     for (PetscInt r = 0; flg && r < nr; r++)
37069d0448ceSStefano Zampini       for (PetscInt k = r + 1; flg && k < nr; k++)
37079d0448ceSStefano Zampini         if ((m[k] <= m[r] && m[r] <= M[k]) || (m[k] <= M[r] && M[r] <= M[k])) flg = PETSC_FALSE;
37089d0448ceSStefano Zampini     PetscCall(PetscFree2(m, M));
37099d0448ceSStefano Zampini     PetscCall(PetscFree2(rows, cols));
371003e5aca4SStefano Zampini     if (!flg) {
371103e5aca4SStefano Zampini       PetscCall(PetscInfo(A, "MAT_FACTOR_CHOLESKY not supported for intersecting row maps. Use MAT_FACTOR_LU.\n"));
371203e5aca4SStefano Zampini       PetscFunctionReturn(PETSC_SUCCESS);
371303e5aca4SStefano Zampini     }
37149d0448ceSStefano Zampini   }
37159d0448ceSStefano Zampini 
37169d0448ceSStefano Zampini   for (PetscInt r = 0; r < nr; r++) {
37179d0448ceSStefano Zampini     for (PetscInt c = 0; c < nc; c++) {
37189d0448ceSStefano Zampini       Mat       sub = mats[r][c];
37194b9405b2SPierre Jolivet       PetscBool isSeqAIJ, isMPIAIJ, isSeqBAIJ, isMPIBAIJ, isSeqSBAIJ, isMPISBAIJ, isTrans, isDiag, isDense;
37209d0448ceSStefano Zampini 
37219d0448ceSStefano Zampini       if (!sub || (ftype == MAT_FACTOR_CHOLESKY && c < r)) continue;
37225d955bbbSStefano Zampini       PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATTRANSPOSEVIRTUAL, &isTrans));
37235d955bbbSStefano Zampini       if (isTrans) PetscCall(MatTransposeGetMat(sub, &sub));
37245d955bbbSStefano Zampini       else {
37255d955bbbSStefano Zampini         PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATHERMITIANTRANSPOSEVIRTUAL, &isTrans));
37265d955bbbSStefano Zampini         if (isTrans) PetscCall(MatHermitianTransposeGetMat(sub, &sub));
37275d955bbbSStefano Zampini       }
37289d0448ceSStefano Zampini       PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATSEQAIJ, &isSeqAIJ));
37299d0448ceSStefano Zampini       PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATMPIAIJ, &isMPIAIJ));
37309d0448ceSStefano Zampini       PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATSEQBAIJ, &isSeqBAIJ));
37319d0448ceSStefano Zampini       PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATMPIBAIJ, &isMPIBAIJ));
37329d0448ceSStefano Zampini       PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATSEQSBAIJ, &isSeqSBAIJ));
37339d0448ceSStefano Zampini       PetscCall(PetscObjectBaseTypeCompare((PetscObject)sub, MATMPISBAIJ, &isMPISBAIJ));
3734c3e1b152SPierre Jolivet       PetscCall(PetscObjectTypeCompare((PetscObject)sub, MATDIAGONAL, &isDiag));
37354b9405b2SPierre Jolivet       PetscCall(PetscObjectTypeCompareAny((PetscObject)sub, &isDense, MATSEQDENSE, MATMPIDENSE, NULL));
37369d0448ceSStefano Zampini       if (ftype == MAT_FACTOR_CHOLESKY) {
3737dcab004fSPierre Jolivet         if (r == c) {
37384b9405b2SPierre Jolivet           if (!isSeqAIJ && !isMPIAIJ && !isSeqBAIJ && !isMPIBAIJ && !isSeqSBAIJ && !isMPISBAIJ && !isDiag && !isDense) {
373940afc089SBarry Smith             PetscCall(PetscInfo(sub, "MAT_FACTOR_CHOLESKY not supported for diagonal block of type %s.\n", ((PetscObject)sub)->type_name));
374003e5aca4SStefano Zampini             flg = PETSC_FALSE;
3741dcab004fSPierre Jolivet           }
37424b9405b2SPierre Jolivet         } else if (!isSeqAIJ && !isMPIAIJ && !isSeqBAIJ && !isMPIBAIJ && !isDiag && !isDense) {
374340afc089SBarry Smith           PetscCall(PetscInfo(sub, "MAT_FACTOR_CHOLESKY not supported for off-diagonal block of type %s.\n", ((PetscObject)sub)->type_name));
374403e5aca4SStefano Zampini           flg = PETSC_FALSE;
374503e5aca4SStefano Zampini         }
37464b9405b2SPierre Jolivet       } else if (!isSeqAIJ && !isMPIAIJ && !isSeqBAIJ && !isMPIBAIJ && !isDiag && !isDense) {
37479afb9c56SPierre Jolivet         PetscCall(PetscInfo(sub, "MAT_FACTOR_LU not supported for block of type %s.\n", ((PetscObject)sub)->type_name));
374803e5aca4SStefano Zampini         flg = PETSC_FALSE;
37499d0448ceSStefano Zampini       }
37509d0448ceSStefano Zampini     }
375103e5aca4SStefano Zampini   }
375203e5aca4SStefano Zampini   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
37539d0448ceSStefano Zampini 
37549d0448ceSStefano Zampini   /* Create the factorization matrix */
37559d0448ceSStefano Zampini   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
37569d0448ceSStefano Zampini   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
37579d0448ceSStefano Zampini   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &((PetscObject)B)->type_name));
37589d0448ceSStefano Zampini   PetscCall(MatSetUp(B));
37599d0448ceSStefano Zampini 
37609d0448ceSStefano Zampini   PetscCall(PetscNew(&mumps));
37619d0448ceSStefano Zampini 
37629d0448ceSStefano Zampini   B->ops->view    = MatView_MUMPS;
37639d0448ceSStefano Zampini   B->ops->getinfo = MatGetInfo_MUMPS;
37649d0448ceSStefano Zampini 
37659d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mumps));
37669d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MUMPS));
37679d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorCreateSchurComplement_C", MatFactorCreateSchurComplement_MUMPS));
37689d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetIcntl_C", MatMumpsSetIcntl_MUMPS));
37699d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetIcntl_C", MatMumpsGetIcntl_MUMPS));
37709d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsSetCntl_C", MatMumpsSetCntl_MUMPS));
37719d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetCntl_C", MatMumpsGetCntl_MUMPS));
37729d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfo_C", MatMumpsGetInfo_MUMPS));
37739d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInfog_C", MatMumpsGetInfog_MUMPS));
37749d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfo_C", MatMumpsGetRinfo_MUMPS));
37759d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetRinfog_C", MatMumpsGetRinfog_MUMPS));
37769d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetNullPivots_C", MatMumpsGetNullPivots_MUMPS));
37779d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverse_C", MatMumpsGetInverse_MUMPS));
37789d0448ceSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMumpsGetInverseTranspose_C", MatMumpsGetInverseTranspose_MUMPS));
37799d0448ceSStefano Zampini 
37809d0448ceSStefano Zampini   if (ftype == MAT_FACTOR_LU) {
37819d0448ceSStefano Zampini     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
37829d0448ceSStefano Zampini     B->factortype            = MAT_FACTOR_LU;
37839d0448ceSStefano Zampini     mumps->sym               = 0;
37849d0448ceSStefano Zampini   } else {
37859d0448ceSStefano Zampini     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
37869d0448ceSStefano Zampini     B->factortype                  = MAT_FACTOR_CHOLESKY;
37879d0448ceSStefano Zampini #if defined(PETSC_USE_COMPLEX)
37889d0448ceSStefano Zampini     mumps->sym = 2;
37899d0448ceSStefano Zampini #else
37909d0448ceSStefano Zampini     if (A->spd == PETSC_BOOL3_TRUE) mumps->sym = 1;
37919d0448ceSStefano Zampini     else mumps->sym = 2;
37929d0448ceSStefano Zampini #endif
37939d0448ceSStefano Zampini   }
37949d0448ceSStefano Zampini   mumps->ConvertToTriples = MatConvertToTriples_nest_xaij;
37959d0448ceSStefano Zampini   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[ftype]));
37969d0448ceSStefano Zampini 
37979d0448ceSStefano Zampini   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
37989d0448ceSStefano Zampini   if (size == 1) {
37999d0448ceSStefano Zampini     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
38009d0448ceSStefano Zampini     B->canuseordering = PETSC_TRUE;
38019d0448ceSStefano Zampini   }
38029d0448ceSStefano Zampini 
38039d0448ceSStefano Zampini   /* set solvertype */
38049d0448ceSStefano Zampini   PetscCall(PetscFree(B->solvertype));
38059d0448ceSStefano Zampini   PetscCall(PetscStrallocpy(MATSOLVERMUMPS, &B->solvertype));
38069d0448ceSStefano Zampini   B->ops->destroy = MatDestroy_MUMPS;
38079d0448ceSStefano Zampini   B->data         = (void *)mumps;
38089d0448ceSStefano Zampini 
38099d0448ceSStefano Zampini   *F               = B;
38109d0448ceSStefano Zampini   mumps->id.job    = JOB_NULL;
38119d0448ceSStefano Zampini   mumps->ICNTL_pre = NULL;
38129d0448ceSStefano Zampini   mumps->CNTL_pre  = NULL;
38139d0448ceSStefano Zampini   mumps->matstruc  = DIFFERENT_NONZERO_PATTERN;
38149d0448ceSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
38159d0448ceSStefano Zampini }
38169d0448ceSStefano Zampini 
3817d1f0640dSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
3818d71ae5a4SJacob Faibussowitsch {
381942c9c57cSBarry Smith   PetscFunctionBegin;
38209566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
38219566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
38229566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIBAIJ, MAT_FACTOR_LU, MatGetFactor_baij_mumps));
38239566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_baij_mumps));
38249566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPISBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_sbaij_mumps));
38259566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
38269566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
38279566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQBAIJ, MAT_FACTOR_LU, MatGetFactor_baij_mumps));
38289566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_baij_mumps));
38299566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_sbaij_mumps));
38309566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQSELL, MAT_FACTOR_LU, MatGetFactor_sell_mumps));
3831c3e1b152SPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATDIAGONAL, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
3832c3e1b152SPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATDIAGONAL, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
38334b9405b2SPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQDENSE, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
38344b9405b2SPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATSEQDENSE, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
38354b9405b2SPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIDENSE, MAT_FACTOR_LU, MatGetFactor_aij_mumps));
38364b9405b2SPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATMPIDENSE, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mumps));
38379d0448ceSStefano Zampini   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATNEST, MAT_FACTOR_LU, MatGetFactor_nest_mumps));
38389d0448ceSStefano Zampini   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS, MATNEST, MAT_FACTOR_CHOLESKY, MatGetFactor_nest_mumps));
38393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
384042c9c57cSBarry Smith }
3841