xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
567602552SJunchao Zhang #include <petscpkg_version.h>
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
87ee00b23SStefano Zampini #include <../src/mat/impls/sell/mpi/mpisell.h>
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
122907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
132907cef9SHong Zhang #include <cmumps_c.h>
142907cef9SHong Zhang #else
15c6db04a5SJed Brown #include <zmumps_c.h>
162907cef9SHong Zhang #endif
172907cef9SHong Zhang #else
182907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
192907cef9SHong Zhang #include <smumps_c.h>
20397b6df1SKris Buschelman #else
21c6db04a5SJed Brown #include <dmumps_c.h>
22397b6df1SKris Buschelman #endif
232907cef9SHong Zhang #endif
24397b6df1SKris Buschelman EXTERN_C_END
25397b6df1SKris Buschelman #define JOB_INIT -1
263d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
273d472b54SHong Zhang #define JOB_FACTNUMERIC 2
283d472b54SHong Zhang #define JOB_SOLVE 3
29397b6df1SKris Buschelman #define JOB_END -2
303d472b54SHong Zhang 
312907cef9SHong Zhang /* calls to MUMPS */
322907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
332907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
343ab56b82SJunchao Zhang #define MUMPS_c cmumps_c
352907cef9SHong Zhang #else
363ab56b82SJunchao Zhang #define MUMPS_c zmumps_c
372907cef9SHong Zhang #endif
382907cef9SHong Zhang #else
392907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
403ab56b82SJunchao Zhang #define MUMPS_c smumps_c
412907cef9SHong Zhang #else
423ab56b82SJunchao Zhang #define MUMPS_c dmumps_c
432907cef9SHong Zhang #endif
442907cef9SHong Zhang #endif
452907cef9SHong Zhang 
46a6053eceSJunchao Zhang /* MUMPS uses MUMPS_INT for nonzero indices such as irn/jcn, irn_loc/jcn_loc and uses int64_t for
47a6053eceSJunchao Zhang    number of nonzeros such as nnz, nnz_loc. We typedef MUMPS_INT to PetscMUMPSInt to follow the
48a6053eceSJunchao Zhang    naming convention in PetscMPIInt, PetscBLASInt etc.
49a6053eceSJunchao Zhang */
50a6053eceSJunchao Zhang typedef MUMPS_INT PetscMUMPSInt;
51a6053eceSJunchao Zhang 
5267602552SJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_GE(5,3,0)
5367602552SJunchao 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 */
54a6053eceSJunchao Zhang     #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
5567602552SJunchao Zhang   #endif
56a6053eceSJunchao Zhang #else
5767602552SJunchao Zhang   #if defined(INTSIZE64) /* INTSIZE64 is a command line macro one used to build MUMPS in full 64-bit mode */
5867602552SJunchao Zhang     #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
5967602552SJunchao Zhang   #endif
6067602552SJunchao Zhang #endif
6167602552SJunchao Zhang 
62a6053eceSJunchao Zhang #define MPIU_MUMPSINT             MPI_INT
63a6053eceSJunchao Zhang #define PETSC_MUMPS_INT_MAX       2147483647
64a6053eceSJunchao Zhang #define PETSC_MUMPS_INT_MIN       -2147483648
65a6053eceSJunchao Zhang 
66a6053eceSJunchao Zhang /* Cast PetscInt to PetscMUMPSInt. Usually there is no overflow since <a> is row/col indices or some small integers*/
679fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscMUMPSIntCast(PetscInt a,PetscMUMPSInt *b)
68a6053eceSJunchao Zhang {
69a6053eceSJunchao Zhang   PetscFunctionBegin;
70ece88022SPierre Jolivet #if PetscDefined(USE_64BIT_INDICES)
712c71b3e2SJacob Faibussowitsch   PetscAssert(a <= PETSC_MUMPS_INT_MAX && a >= PETSC_MUMPS_INT_MIN,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PetscInt too long for PetscMUMPSInt");
72ece88022SPierre Jolivet #endif
73a6053eceSJunchao Zhang   *b = (PetscMUMPSInt)(a);
74a6053eceSJunchao Zhang   PetscFunctionReturn(0);
75a6053eceSJunchao Zhang }
76a6053eceSJunchao Zhang 
77a6053eceSJunchao Zhang /* Put these utility routines here since they are only used in this file */
789fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscOptionsMUMPSInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscMUMPSInt currentvalue,PetscMUMPSInt *value,PetscBool *set,PetscMUMPSInt lb,PetscMUMPSInt ub)
79a6053eceSJunchao Zhang {
80a6053eceSJunchao Zhang   PetscInt       myval;
81a6053eceSJunchao Zhang   PetscBool      myset;
82a6053eceSJunchao Zhang   PetscFunctionBegin;
83a6053eceSJunchao Zhang   /* PetscInt's size should be always >= PetscMUMPSInt's. It is safe to call PetscOptionsInt_Private to read a PetscMUMPSInt */
849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt_Private(PetscOptionsObject,opt,text,man,(PetscInt)currentvalue,&myval,&myset,lb,ub));
859566063dSJacob Faibussowitsch   if (myset) PetscCall(PetscMUMPSIntCast(myval,value));
86a6053eceSJunchao Zhang   if (set) *set = myset;
87a6053eceSJunchao Zhang   PetscFunctionReturn(0);
88a6053eceSJunchao Zhang }
89a6053eceSJunchao 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)
90a6053eceSJunchao Zhang 
91217d3b1eSJunchao 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 */
923ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
933ab56b82SJunchao Zhang #define PetscMUMPS_c(mumps) \
943ab56b82SJunchao Zhang   do { \
953ab56b82SJunchao Zhang     if (mumps->use_petsc_omp_support) { \
963ab56b82SJunchao Zhang       if (mumps->is_omp_master) { \
979566063dSJacob Faibussowitsch         PetscCall(PetscOmpCtrlOmpRegionOnMasterBegin(mumps->omp_ctrl)); \
983ab56b82SJunchao Zhang         MUMPS_c(&mumps->id); \
999566063dSJacob Faibussowitsch         PetscCall(PetscOmpCtrlOmpRegionOnMasterEnd(mumps->omp_ctrl)); \
1003ab56b82SJunchao Zhang       } \
1019566063dSJacob Faibussowitsch       PetscCall(PetscOmpCtrlBarrier(mumps->omp_ctrl)); \
102c3714a1dSJunchao Zhang       /* Global info is same on all processes so we Bcast it within omp_comm. Local info is specific      \
103c3714a1dSJunchao Zhang          to processes, so we only Bcast info[1], an error code and leave others (since they do not have   \
104c3714a1dSJunchao Zhang          an easy translation between omp_comm and petsc_comm). See MUMPS-5.1.2 manual p82.                   \
105c3714a1dSJunchao Zhang          omp_comm is a small shared memory communicator, hence doing multiple Bcast as shown below is OK. \
106c3714a1dSJunchao Zhang       */ \
1079566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(mumps->id.infog, 40,MPIU_MUMPSINT, 0,mumps->omp_comm));\
1089566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(mumps->id.rinfog,20,MPIU_REAL,     0,mumps->omp_comm));\
1099566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(mumps->id.info,  1, MPIU_MUMPSINT, 0,mumps->omp_comm));\
1103ab56b82SJunchao Zhang     } else { \
1113ab56b82SJunchao Zhang       MUMPS_c(&mumps->id); \
1123ab56b82SJunchao Zhang     } \
1133ab56b82SJunchao Zhang   } while (0)
1143ab56b82SJunchao Zhang #else
1153ab56b82SJunchao Zhang #define PetscMUMPS_c(mumps) \
1163ab56b82SJunchao Zhang   do { MUMPS_c(&mumps->id); } while (0)
1173ab56b82SJunchao Zhang #endif
1183ab56b82SJunchao Zhang 
119940cd9d6SSatish Balay /* declare MumpsScalar */
120940cd9d6SSatish Balay #if defined(PETSC_USE_COMPLEX)
121940cd9d6SSatish Balay #if defined(PETSC_USE_REAL_SINGLE)
122940cd9d6SSatish Balay #define MumpsScalar mumps_complex
123940cd9d6SSatish Balay #else
124940cd9d6SSatish Balay #define MumpsScalar mumps_double_complex
125940cd9d6SSatish Balay #endif
126940cd9d6SSatish Balay #else
127940cd9d6SSatish Balay #define MumpsScalar PetscScalar
128940cd9d6SSatish Balay #endif
1293d472b54SHong Zhang 
130397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
131397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
132397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
133397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
134a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
135397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
136adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
137397b6df1SKris Buschelman 
138a6053eceSJunchao Zhang typedef struct Mat_MUMPS Mat_MUMPS;
139a6053eceSJunchao Zhang struct Mat_MUMPS {
140397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
1412907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1422907cef9SHong Zhang   CMUMPS_STRUC_C id;
1432907cef9SHong Zhang #else
144397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
1452907cef9SHong Zhang #endif
1462907cef9SHong Zhang #else
1472907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1482907cef9SHong Zhang   SMUMPS_STRUC_C id;
149397b6df1SKris Buschelman #else
150397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
151397b6df1SKris Buschelman #endif
1522907cef9SHong Zhang #endif
1532907cef9SHong Zhang 
154397b6df1SKris Buschelman   MatStructure   matstruc;
1552d4298aeSJunchao Zhang   PetscMPIInt    myid,petsc_size;
156a6053eceSJunchao Zhang   PetscMUMPSInt  *irn,*jcn;             /* the (i,j,v) triplets passed to mumps. */
157a6053eceSJunchao 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. */
158a6053eceSJunchao Zhang   PetscInt64     nnz;                   /* number of nonzeros. The type is called selective 64-bit in mumps */
159a6053eceSJunchao Zhang   PetscMUMPSInt  sym;
1602d4298aeSJunchao Zhang   MPI_Comm       mumps_comm;
161a6053eceSJunchao Zhang   PetscMUMPSInt  ICNTL9_pre;            /* check if ICNTL(9) is changed from previous MatSolve */
162801fbe65SHong Zhang   VecScatter     scat_rhs, scat_sol;    /* used by MatSolve() */
16325aac85cSJunchao Zhang   PetscMUMPSInt  ICNTL20;               /* use centralized (0) or distributed (10) dense RHS */
16467602552SJunchao Zhang   PetscMUMPSInt  lrhs_loc,nloc_rhs,*irhs_loc;
16567602552SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
16667602552SJunchao Zhang   PetscInt       *rhs_nrow,max_nrhs;
16767602552SJunchao Zhang   PetscMPIInt    *rhs_recvcounts,*rhs_disps;
16867602552SJunchao Zhang   PetscScalar    *rhs_loc,*rhs_recvbuf;
16967602552SJunchao Zhang #endif
170801fbe65SHong Zhang   Vec            b_seq,x_seq;
171a6053eceSJunchao Zhang   PetscInt       ninfo,*info;           /* which INFO to display */
172b5fa320bSStefano Zampini   PetscInt       sizeredrhs;
17359ac8732SStefano Zampini   PetscScalar    *schur_sol;
17459ac8732SStefano Zampini   PetscInt       schur_sizesol;
175a6053eceSJunchao Zhang   PetscMUMPSInt  *ia_alloc,*ja_alloc;   /* work arrays used for the CSR struct for sparse rhs */
176a6053eceSJunchao Zhang   PetscInt64     cur_ilen,cur_jlen;     /* current len of ia_alloc[], ja_alloc[] */
177a6053eceSJunchao Zhang   PetscErrorCode (*ConvertToTriples)(Mat,PetscInt,MatReuse,Mat_MUMPS*);
1782205254eSKarl Rupp 
179a6053eceSJunchao Zhang   /* stuff used by petsc/mumps OpenMP support*/
1803ab56b82SJunchao Zhang   PetscBool      use_petsc_omp_support;
1813ab56b82SJunchao Zhang   PetscOmpCtrl   omp_ctrl;              /* an OpenMP controler that blocked processes will release their CPU (MPI_Barrier does not have this guarantee) */
1823ab56b82SJunchao Zhang   MPI_Comm       petsc_comm,omp_comm;   /* petsc_comm is petsc matrix's comm */
183a6053eceSJunchao Zhang   PetscInt64     *recvcount;            /* a collection of nnz on omp_master */
184a6053eceSJunchao Zhang   PetscMPIInt    tag,omp_comm_size;
1853ab56b82SJunchao Zhang   PetscBool      is_omp_master;         /* is this rank the master of omp_comm */
186a6053eceSJunchao Zhang   MPI_Request    *reqs;
187a6053eceSJunchao Zhang };
1883ab56b82SJunchao Zhang 
189a6053eceSJunchao Zhang /* Cast a 1-based CSR represented by (nrow, ia, ja) of type PetscInt to a CSR of type PetscMUMPSInt.
190a6053eceSJunchao Zhang    Here, nrow is number of rows, ia[] is row pointer and ja[] is column indices.
191a6053eceSJunchao Zhang  */
192a6053eceSJunchao Zhang static PetscErrorCode PetscMUMPSIntCSRCast(Mat_MUMPS *mumps,PetscInt nrow,PetscInt *ia,PetscInt *ja,PetscMUMPSInt **ia_mumps,PetscMUMPSInt **ja_mumps,PetscMUMPSInt *nnz_mumps)
193a6053eceSJunchao Zhang {
194a6053eceSJunchao Zhang   PetscInt       nnz=ia[nrow]-1; /* mumps uses 1-based indices. Uses PetscInt instead of PetscInt64 since mumps only uses PetscMUMPSInt for rhs */
195f0c56d0fSKris Buschelman 
196a6053eceSJunchao Zhang   PetscFunctionBegin;
197a6053eceSJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
198a6053eceSJunchao Zhang   {
199a6053eceSJunchao Zhang     PetscInt i;
200a6053eceSJunchao Zhang     if (nrow+1 > mumps->cur_ilen) { /* realloc ia_alloc/ja_alloc to fit ia/ja */
2019566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->ia_alloc));
2029566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrow+1,&mumps->ia_alloc));
203a6053eceSJunchao Zhang       mumps->cur_ilen = nrow+1;
204a6053eceSJunchao Zhang     }
205a6053eceSJunchao Zhang     if (nnz > mumps->cur_jlen) {
2069566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->ja_alloc));
2079566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nnz,&mumps->ja_alloc));
208a6053eceSJunchao Zhang       mumps->cur_jlen = nnz;
209a6053eceSJunchao Zhang     }
2109566063dSJacob Faibussowitsch     for (i=0; i<nrow+1; i++) PetscCall(PetscMUMPSIntCast(ia[i],&(mumps->ia_alloc[i])));
2119566063dSJacob Faibussowitsch     for (i=0; i<nnz; i++)    PetscCall(PetscMUMPSIntCast(ja[i],&(mumps->ja_alloc[i])));
212a6053eceSJunchao Zhang     *ia_mumps = mumps->ia_alloc;
213a6053eceSJunchao Zhang     *ja_mumps = mumps->ja_alloc;
214a6053eceSJunchao Zhang   }
215a6053eceSJunchao Zhang #else
216a6053eceSJunchao Zhang   *ia_mumps = ia;
217a6053eceSJunchao Zhang   *ja_mumps = ja;
218a6053eceSJunchao Zhang #endif
2199566063dSJacob Faibussowitsch   PetscCall(PetscMUMPSIntCast(nnz,nnz_mumps));
220a6053eceSJunchao Zhang   PetscFunctionReturn(0);
221a6053eceSJunchao Zhang }
222b24902e0SBarry Smith 
22359ac8732SStefano Zampini static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
224b5fa320bSStefano Zampini {
225b5fa320bSStefano Zampini   PetscFunctionBegin;
2269566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.listvar_schur));
2279566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.redrhs));
2289566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->schur_sol));
22959ac8732SStefano Zampini   mumps->id.size_schur = 0;
230b3cb21ddSStefano Zampini   mumps->id.schur_lld  = 0;
23159ac8732SStefano Zampini   mumps->id.ICNTL(19)  = 0;
23259ac8732SStefano Zampini   PetscFunctionReturn(0);
23359ac8732SStefano Zampini }
23459ac8732SStefano Zampini 
235b3cb21ddSStefano Zampini /* solve with rhs in mumps->id.redrhs and return in the same location */
236b3cb21ddSStefano Zampini static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
23759ac8732SStefano Zampini {
238b3cb21ddSStefano Zampini   Mat_MUMPS            *mumps=(Mat_MUMPS*)F->data;
239b3cb21ddSStefano Zampini   Mat                  S,B,X;
240b3cb21ddSStefano Zampini   MatFactorSchurStatus schurstatus;
241b3cb21ddSStefano Zampini   PetscInt             sizesol;
24259ac8732SStefano Zampini 
24359ac8732SStefano Zampini   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCall(MatFactorFactorizeSchurComplement(F));
2459566063dSJacob Faibussowitsch   PetscCall(MatFactorGetSchurComplement(F,&S,&schurstatus));
2469566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B));
2479566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,((PetscObject)S)->type_name));
248a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2499566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(B,S->boundtocpu));
250a3d589ffSStefano Zampini #endif
251b3cb21ddSStefano Zampini   switch (schurstatus) {
252b3cb21ddSStefano Zampini   case MAT_FACTOR_SCHUR_FACTORED:
2539566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X));
2549566063dSJacob Faibussowitsch     PetscCall(MatSetType(X,((PetscObject)S)->type_name));
255a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2569566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(X,S->boundtocpu));
257a3d589ffSStefano Zampini #endif
258b3cb21ddSStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
2599566063dSJacob Faibussowitsch       PetscCall(MatMatSolveTranspose(S,B,X));
26059ac8732SStefano Zampini     } else {
2619566063dSJacob Faibussowitsch       PetscCall(MatMatSolve(S,B,X));
26259ac8732SStefano Zampini     }
263b3cb21ddSStefano Zampini     break;
264b3cb21ddSStefano Zampini   case MAT_FACTOR_SCHUR_INVERTED:
265b3cb21ddSStefano Zampini     sizesol = mumps->id.nrhs*mumps->id.size_schur;
26659ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
2679566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->schur_sol));
2689566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sizesol,&mumps->schur_sol));
26959ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
270b5fa320bSStefano Zampini     }
2719566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X));
2729566063dSJacob Faibussowitsch     PetscCall(MatSetType(X,((PetscObject)S)->type_name));
273a3d589ffSStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2749566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(X,S->boundtocpu));
275a3d589ffSStefano Zampini #endif
2769566063dSJacob Faibussowitsch     PetscCall(MatProductCreateWithMat(S,B,NULL,X));
27759ac8732SStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
2789566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(X,MATPRODUCT_AtB));
279b5fa320bSStefano Zampini     } else {
2809566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(X,MATPRODUCT_AB));
281b5fa320bSStefano Zampini     }
2829566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(X));
2839566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(X));
2849566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(X));
2854417c5e8SHong Zhang 
2869566063dSJacob Faibussowitsch     PetscCall(MatCopy(X,B,SAME_NONZERO_PATTERN));
287b3cb21ddSStefano Zampini     break;
288b3cb21ddSStefano Zampini   default:
28998921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %d",F->schur_status);
29059ac8732SStefano Zampini   }
2919566063dSJacob Faibussowitsch   PetscCall(MatFactorRestoreSchurComplement(F,&S,schurstatus));
2929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
294b5fa320bSStefano Zampini   PetscFunctionReturn(0);
295b5fa320bSStefano Zampini }
296b5fa320bSStefano Zampini 
297b3cb21ddSStefano Zampini static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
298b5fa320bSStefano Zampini {
299b3cb21ddSStefano Zampini   Mat_MUMPS     *mumps=(Mat_MUMPS*)F->data;
300b5fa320bSStefano Zampini 
301b5fa320bSStefano Zampini   PetscFunctionBegin;
302b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
303b5fa320bSStefano Zampini     PetscFunctionReturn(0);
304b5fa320bSStefano Zampini   }
305b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
306b5fa320bSStefano Zampini     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
307b5fa320bSStefano Zampini     /* allocate MUMPS internal array to store reduced right-hand sides */
308b5fa320bSStefano Zampini     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
3099566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->id.redrhs));
310b5fa320bSStefano Zampini       mumps->id.lredrhs = mumps->id.size_schur;
3119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs));
312b5fa320bSStefano Zampini       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
313b5fa320bSStefano Zampini     }
314b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 1; /* condensation phase */
315b5fa320bSStefano Zampini   } else { /* prepare for the expansion step */
316b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
3179566063dSJacob Faibussowitsch     PetscCall(MatMumpsSolveSchur_Private(F));
318b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
3193ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
320*08401ef6SPierre Jolivet     PetscCheck(mumps->id.INFOG(1) >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d",mumps->id.INFOG(1));
321b5fa320bSStefano Zampini     /* restore defaults */
322b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
323d3d598ffSStefano Zampini     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
324d3d598ffSStefano Zampini     if (mumps->id.nrhs > 1) {
3259566063dSJacob Faibussowitsch       PetscCall(PetscFree(mumps->id.redrhs));
326d3d598ffSStefano Zampini       mumps->id.lredrhs = 0;
327d3d598ffSStefano Zampini       mumps->sizeredrhs = 0;
328d3d598ffSStefano Zampini     }
329b5fa320bSStefano Zampini   }
330b5fa320bSStefano Zampini   PetscFunctionReturn(0);
331b5fa320bSStefano Zampini }
332b5fa320bSStefano Zampini 
333397b6df1SKris Buschelman /*
334d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
335d341cd04SHong Zhang 
336397b6df1SKris Buschelman   input:
33775480915SPierre Jolivet     A       - matrix in aij,baij or sbaij format
338397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
339bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
340bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
341397b6df1SKris Buschelman   output:
342397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
343397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
344eb9baa12SBarry Smith 
345eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
3467ee00b23SStefano Zampini   freed with PetscFree(mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
347eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
348eb9baa12SBarry Smith 
349397b6df1SKris Buschelman  */
35016ebf90aSShri Abhyankar 
351a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
352b24902e0SBarry Smith {
353a3d589ffSStefano Zampini   const PetscScalar *av;
354185f6596SHong Zhang   const PetscInt    *ai,*aj,*ajj,M=A->rmap->n;
355a6053eceSJunchao Zhang   PetscInt64        nz,rnz,i,j,k;
356a6053eceSJunchao Zhang   PetscMUMPSInt     *row,*col;
35716ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
358397b6df1SKris Buschelman 
359397b6df1SKris Buschelman   PetscFunctionBegin;
3609566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(A,&av));
361a6053eceSJunchao Zhang   mumps->val = (PetscScalar*)av;
362bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
3632205254eSKarl Rupp     nz   = aa->nz;
3642205254eSKarl Rupp     ai   = aa->i;
3652205254eSKarl Rupp     aj   = aa->j;
3669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
367a6053eceSJunchao Zhang     for (i=k=0; i<M; i++) {
36816ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
36967877ebaSShri Abhyankar       ajj = aj + ai[i];
37067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
3719566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(i+shift,&row[k]));
3729566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(ajj[j] + shift,&col[k]));
373a6053eceSJunchao Zhang         k++;
37416ebf90aSShri Abhyankar       }
37516ebf90aSShri Abhyankar     }
376a6053eceSJunchao Zhang     mumps->irn = row;
377a6053eceSJunchao Zhang     mumps->jcn = col;
378a6053eceSJunchao Zhang     mumps->nnz = nz;
37916ebf90aSShri Abhyankar   }
3809566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(A,&av));
38116ebf90aSShri Abhyankar   PetscFunctionReturn(0);
38216ebf90aSShri Abhyankar }
383397b6df1SKris Buschelman 
384a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
3857ee00b23SStefano Zampini {
386a6053eceSJunchao Zhang   PetscInt64     nz,i,j,k,r;
3877ee00b23SStefano Zampini   Mat_SeqSELL    *a=(Mat_SeqSELL*)A->data;
388a6053eceSJunchao Zhang   PetscMUMPSInt  *row,*col;
3897ee00b23SStefano Zampini 
3907ee00b23SStefano Zampini   PetscFunctionBegin;
391a6053eceSJunchao Zhang   mumps->val = a->val;
3927ee00b23SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
3937ee00b23SStefano Zampini     nz   = a->sliidx[a->totalslices];
3949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
395a6053eceSJunchao Zhang     for (i=k=0; i<a->totalslices; i++) {
396a6053eceSJunchao Zhang       for (j=a->sliidx[i],r=0; j<a->sliidx[i+1]; j++,r=((r+1)&0x07)) {
3979566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(8*i+r+shift,&row[k++]));
3987ee00b23SStefano Zampini       }
3997ee00b23SStefano Zampini     }
4009566063dSJacob Faibussowitsch     for (i=0;i<nz;i++) PetscCall(PetscMUMPSIntCast(a->colidx[i]+shift,&col[i]));
401a6053eceSJunchao Zhang     mumps->irn = row;
402a6053eceSJunchao Zhang     mumps->jcn = col;
403a6053eceSJunchao Zhang     mumps->nnz = nz;
4047ee00b23SStefano Zampini   }
4057ee00b23SStefano Zampini   PetscFunctionReturn(0);
4067ee00b23SStefano Zampini }
4077ee00b23SStefano Zampini 
408a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
40967877ebaSShri Abhyankar {
41067877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
41133d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
412a6053eceSJunchao Zhang   PetscInt64     M,nz,idx=0,rnz,i,j,k,m;
413a6053eceSJunchao Zhang   PetscInt       bs;
414a6053eceSJunchao Zhang   PetscMUMPSInt  *row,*col;
41567877ebaSShri Abhyankar 
41667877ebaSShri Abhyankar   PetscFunctionBegin;
4179566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
41833d57670SJed Brown   M          = A->rmap->N/bs;
419a6053eceSJunchao Zhang   mumps->val = aa->a;
420bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
421cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
42267877ebaSShri Abhyankar     nz   = bs2*aa->nz;
4239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
42467877ebaSShri Abhyankar     for (i=0; i<M; i++) {
42567877ebaSShri Abhyankar       ajj = aj + ai[i];
42667877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
42767877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
42867877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
42967877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
4309566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i*bs + m + shift,&row[idx]));
4319566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(bs*ajj[k] + j + shift,&col[idx]));
432a6053eceSJunchao Zhang             idx++;
43367877ebaSShri Abhyankar           }
43467877ebaSShri Abhyankar         }
43567877ebaSShri Abhyankar       }
43667877ebaSShri Abhyankar     }
437a6053eceSJunchao Zhang     mumps->irn = row;
438a6053eceSJunchao Zhang     mumps->jcn = col;
439a6053eceSJunchao Zhang     mumps->nnz = nz;
44067877ebaSShri Abhyankar   }
44167877ebaSShri Abhyankar   PetscFunctionReturn(0);
44267877ebaSShri Abhyankar }
44367877ebaSShri Abhyankar 
444a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
44516ebf90aSShri Abhyankar {
44675480915SPierre Jolivet   const PetscInt *ai, *aj,*ajj;
447a6053eceSJunchao Zhang   PetscInt        bs;
448a6053eceSJunchao Zhang   PetscInt64      nz,rnz,i,j,k,m;
449a6053eceSJunchao Zhang   PetscMUMPSInt   *row,*col;
45075480915SPierre Jolivet   PetscScalar     *val;
45116ebf90aSShri Abhyankar   Mat_SeqSBAIJ    *aa=(Mat_SeqSBAIJ*)A->data;
45275480915SPierre Jolivet   const PetscInt  bs2=aa->bs2,mbs=aa->mbs;
45338548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
45438548759SBarry Smith   PetscBool       hermitian;
45538548759SBarry Smith #endif
45616ebf90aSShri Abhyankar 
45716ebf90aSShri Abhyankar   PetscFunctionBegin;
45838548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
4599566063dSJacob Faibussowitsch   PetscCall(MatGetOption(A,MAT_HERMITIAN,&hermitian));
46028b400f6SJacob Faibussowitsch   PetscCheck(!hermitian,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
46138548759SBarry Smith #endif
4622205254eSKarl Rupp   ai   = aa->i;
4632205254eSKarl Rupp   aj   = aa->j;
4649566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
46575480915SPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
46675480915SPierre Jolivet     nz   = aa->nz;
4679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(bs2*nz,&row,bs2*nz,&col));
468a6053eceSJunchao Zhang     if (bs>1) {
4699566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs2*nz,&mumps->val_alloc));
470a6053eceSJunchao Zhang       mumps->val = mumps->val_alloc;
47175480915SPierre Jolivet     } else {
472a6053eceSJunchao Zhang       mumps->val = aa->a;
47375480915SPierre Jolivet     }
474a6053eceSJunchao Zhang     mumps->irn = row;
475a6053eceSJunchao Zhang     mumps->jcn = col;
476a6053eceSJunchao Zhang   } else {
477a6053eceSJunchao Zhang     if (bs == 1) mumps->val = aa->a;
478a6053eceSJunchao Zhang     row = mumps->irn;
479a6053eceSJunchao Zhang     col = mumps->jcn;
480a6053eceSJunchao Zhang   }
481a6053eceSJunchao Zhang   val = mumps->val;
482185f6596SHong Zhang 
48316ebf90aSShri Abhyankar   nz = 0;
484a81fe166SPierre Jolivet   if (bs>1) {
48575480915SPierre Jolivet     for (i=0; i<mbs; i++) {
48616ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
48767877ebaSShri Abhyankar       ajj = aj + ai[i];
48875480915SPierre Jolivet       for (j=0; j<rnz; j++) {
48975480915SPierre Jolivet         for (k=0; k<bs; k++) {
49075480915SPierre Jolivet           for (m=0; m<bs; m++) {
491ec4f40fdSPierre Jolivet             if (ajj[j]>i || k>=m) {
49275480915SPierre Jolivet               if (reuse == MAT_INITIAL_MATRIX) {
4939566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(i*bs + m + shift,&row[nz]));
4949566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(ajj[j]*bs + k + shift,&col[nz]));
49575480915SPierre Jolivet               }
49675480915SPierre Jolivet               val[nz++] = aa->a[(ai[i]+j)*bs2 + m + k*bs];
49775480915SPierre Jolivet             }
49875480915SPierre Jolivet           }
49975480915SPierre Jolivet         }
50075480915SPierre Jolivet       }
50175480915SPierre Jolivet     }
502a81fe166SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) {
503a81fe166SPierre Jolivet     for (i=0; i<mbs; i++) {
504a81fe166SPierre Jolivet       rnz = ai[i+1] - ai[i];
505a81fe166SPierre Jolivet       ajj = aj + ai[i];
506a81fe166SPierre Jolivet       for (j=0; j<rnz; j++) {
5079566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(i+shift,&row[nz]));
5089566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(ajj[j] + shift,&col[nz]));
509a6053eceSJunchao Zhang         nz++;
510a81fe166SPierre Jolivet       }
511a81fe166SPierre Jolivet     }
512*08401ef6SPierre Jolivet     PetscCheck(nz == aa->nz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Different numbers of nonzeros %" PetscInt64_FMT " != %" PetscInt_FMT,nz,aa->nz);
51375480915SPierre Jolivet   }
514a6053eceSJunchao Zhang   if (reuse == MAT_INITIAL_MATRIX) mumps->nnz = nz;
51516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
51616ebf90aSShri Abhyankar }
51716ebf90aSShri Abhyankar 
518a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
51916ebf90aSShri Abhyankar {
52067877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
521a6053eceSJunchao Zhang   PetscInt64        nz,rnz,i,j;
52267877ebaSShri Abhyankar   const PetscScalar *av,*v1;
52316ebf90aSShri Abhyankar   PetscScalar       *val;
524a6053eceSJunchao Zhang   PetscMUMPSInt     *row,*col;
525829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
52629b521d4Sstefano_zampini   PetscBool         missing;
52738548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
52838548759SBarry Smith   PetscBool         hermitian;
52938548759SBarry Smith #endif
53016ebf90aSShri Abhyankar 
53116ebf90aSShri Abhyankar   PetscFunctionBegin;
53238548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
5339566063dSJacob Faibussowitsch   PetscCall(MatGetOption(A,MAT_HERMITIAN,&hermitian));
53428b400f6SJacob Faibussowitsch   PetscCheck(!hermitian,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
53538548759SBarry Smith #endif
5369566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(A,&av));
537a3d589ffSStefano Zampini   ai    = aa->i; aj = aa->j;
53816ebf90aSShri Abhyankar   adiag = aa->diag;
5399566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal_SeqAIJ(A,&missing,NULL));
540bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
5417ee00b23SStefano Zampini     /* count nz in the upper triangular part of A */
542829b1710SHong Zhang     nz = 0;
54329b521d4Sstefano_zampini     if (missing) {
54429b521d4Sstefano_zampini       for (i=0; i<M; i++) {
54529b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
54629b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
54729b521d4Sstefano_zampini             if (aj[j] < i) continue;
54829b521d4Sstefano_zampini             nz++;
54929b521d4Sstefano_zampini           }
55029b521d4Sstefano_zampini         } else {
55129b521d4Sstefano_zampini           nz += ai[i+1] - adiag[i];
55229b521d4Sstefano_zampini         }
55329b521d4Sstefano_zampini       }
55429b521d4Sstefano_zampini     } else {
555829b1710SHong Zhang       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
55629b521d4Sstefano_zampini     }
5579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
5589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz,&val));
559a6053eceSJunchao Zhang     mumps->nnz = nz;
560a6053eceSJunchao Zhang     mumps->irn = row;
561a6053eceSJunchao Zhang     mumps->jcn = col;
562a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
563185f6596SHong Zhang 
56416ebf90aSShri Abhyankar     nz = 0;
56529b521d4Sstefano_zampini     if (missing) {
56629b521d4Sstefano_zampini       for (i=0; i<M; i++) {
56729b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
56829b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
56929b521d4Sstefano_zampini             if (aj[j] < i) continue;
5709566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i+shift,&row[nz]));
5719566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(aj[j]+shift,&col[nz]));
57229b521d4Sstefano_zampini             val[nz] = av[j];
57329b521d4Sstefano_zampini             nz++;
57429b521d4Sstefano_zampini           }
57529b521d4Sstefano_zampini         } else {
57629b521d4Sstefano_zampini           rnz = ai[i+1] - adiag[i];
57729b521d4Sstefano_zampini           ajj = aj + adiag[i];
57829b521d4Sstefano_zampini           v1  = av + adiag[i];
57929b521d4Sstefano_zampini           for (j=0; j<rnz; j++) {
5809566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(i+shift,&row[nz]));
5819566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(ajj[j] + shift,&col[nz]));
582a6053eceSJunchao Zhang             val[nz++] = v1[j];
58329b521d4Sstefano_zampini           }
58429b521d4Sstefano_zampini         }
58529b521d4Sstefano_zampini       }
58629b521d4Sstefano_zampini     } else {
58716ebf90aSShri Abhyankar       for (i=0; i<M; i++) {
58816ebf90aSShri Abhyankar         rnz = ai[i+1] - adiag[i];
58967877ebaSShri Abhyankar         ajj = aj + adiag[i];
590cf3759fdSShri Abhyankar         v1  = av + adiag[i];
59167877ebaSShri Abhyankar         for (j=0; j<rnz; j++) {
5929566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(i+shift,&row[nz]));
5939566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(ajj[j] + shift,&col[nz]));
594a6053eceSJunchao Zhang           val[nz++] = v1[j];
59516ebf90aSShri Abhyankar         }
59616ebf90aSShri Abhyankar       }
59729b521d4Sstefano_zampini     }
598397b6df1SKris Buschelman   } else {
599a6053eceSJunchao Zhang     nz = 0;
600a6053eceSJunchao Zhang     val = mumps->val;
60129b521d4Sstefano_zampini     if (missing) {
60216ebf90aSShri Abhyankar       for (i=0; i <M; i++) {
60329b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
60429b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
60529b521d4Sstefano_zampini             if (aj[j] < i) continue;
60629b521d4Sstefano_zampini             val[nz++] = av[j];
60729b521d4Sstefano_zampini           }
60829b521d4Sstefano_zampini         } else {
60916ebf90aSShri Abhyankar           rnz = ai[i+1] - adiag[i];
61067877ebaSShri Abhyankar           v1  = av + adiag[i];
61167877ebaSShri Abhyankar           for (j=0; j<rnz; j++) {
61267877ebaSShri Abhyankar             val[nz++] = v1[j];
61316ebf90aSShri Abhyankar           }
61416ebf90aSShri Abhyankar         }
61516ebf90aSShri Abhyankar       }
61629b521d4Sstefano_zampini     } else {
61716ebf90aSShri Abhyankar       for (i=0; i <M; i++) {
61816ebf90aSShri Abhyankar         rnz = ai[i+1] - adiag[i];
61916ebf90aSShri Abhyankar         v1  = av + adiag[i];
62016ebf90aSShri Abhyankar         for (j=0; j<rnz; j++) {
62116ebf90aSShri Abhyankar           val[nz++] = v1[j];
62216ebf90aSShri Abhyankar         }
62316ebf90aSShri Abhyankar       }
62416ebf90aSShri Abhyankar     }
62529b521d4Sstefano_zampini   }
6269566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(A,&av));
62716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
62816ebf90aSShri Abhyankar }
62916ebf90aSShri Abhyankar 
630a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
63116ebf90aSShri Abhyankar {
632a6053eceSJunchao Zhang   const PetscInt    *ai,*aj,*bi,*bj,*garray,*ajj,*bjj;
633a6053eceSJunchao Zhang   PetscInt          bs;
634a6053eceSJunchao Zhang   PetscInt64        rstart,nz,i,j,k,m,jj,irow,countA,countB;
635a6053eceSJunchao Zhang   PetscMUMPSInt     *row,*col;
63616ebf90aSShri Abhyankar   const PetscScalar *av,*bv,*v1,*v2;
63716ebf90aSShri Abhyankar   PetscScalar       *val;
638397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
639397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
640397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
641ec4f40fdSPierre Jolivet   const PetscInt    bs2=aa->bs2,mbs=aa->mbs;
64238548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
64338548759SBarry Smith   PetscBool         hermitian;
64438548759SBarry Smith #endif
64516ebf90aSShri Abhyankar 
64616ebf90aSShri Abhyankar   PetscFunctionBegin;
64738548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
6489566063dSJacob Faibussowitsch   PetscCall(MatGetOption(A,MAT_HERMITIAN,&hermitian));
64928b400f6SJacob Faibussowitsch   PetscCheck(!hermitian,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
65038548759SBarry Smith #endif
6519566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
65238548759SBarry Smith   rstart = A->rmap->rstart;
65338548759SBarry Smith   ai = aa->i;
65438548759SBarry Smith   aj = aa->j;
65538548759SBarry Smith   bi = bb->i;
65638548759SBarry Smith   bj = bb->j;
65738548759SBarry Smith   av = aa->a;
65838548759SBarry Smith   bv = bb->a;
659397b6df1SKris Buschelman 
6602205254eSKarl Rupp   garray = mat->garray;
6612205254eSKarl Rupp 
662bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
663a6053eceSJunchao Zhang     nz   = (aa->nz+bb->nz)*bs2; /* just a conservative estimate */
6649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
6659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz,&val));
666a6053eceSJunchao Zhang     /* can not decide the exact mumps->nnz now because of the SBAIJ */
667a6053eceSJunchao Zhang     mumps->irn = row;
668a6053eceSJunchao Zhang     mumps->jcn = col;
669a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
670397b6df1SKris Buschelman   } else {
671a6053eceSJunchao Zhang     val = mumps->val;
672397b6df1SKris Buschelman   }
673397b6df1SKris Buschelman 
674028e57e8SHong Zhang   jj = 0; irow = rstart;
675ec4f40fdSPierre Jolivet   for (i=0; i<mbs; i++) {
676397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
677397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
678397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
679397b6df1SKris Buschelman     bjj    = bj + bi[i];
680ec4f40fdSPierre Jolivet     v1     = av + ai[i]*bs2;
681ec4f40fdSPierre Jolivet     v2     = bv + bi[i]*bs2;
682397b6df1SKris Buschelman 
683ec4f40fdSPierre Jolivet     if (bs>1) {
684ec4f40fdSPierre Jolivet       /* A-part */
685ec4f40fdSPierre Jolivet       for (j=0; j<countA; j++) {
686ec4f40fdSPierre Jolivet         for (k=0; k<bs; k++) {
687ec4f40fdSPierre Jolivet           for (m=0; m<bs; m++) {
688ec4f40fdSPierre Jolivet             if (rstart + ajj[j]*bs>irow || k>=m) {
689ec4f40fdSPierre Jolivet               if (reuse == MAT_INITIAL_MATRIX) {
6909566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(irow + m + shift,&row[jj]));
6919566063dSJacob Faibussowitsch                 PetscCall(PetscMUMPSIntCast(rstart + ajj[j]*bs + k + shift,&col[jj]));
692ec4f40fdSPierre Jolivet               }
693ec4f40fdSPierre Jolivet               val[jj++] = v1[j*bs2 + m + k*bs];
694ec4f40fdSPierre Jolivet             }
695ec4f40fdSPierre Jolivet           }
696ec4f40fdSPierre Jolivet         }
697ec4f40fdSPierre Jolivet       }
698ec4f40fdSPierre Jolivet 
699ec4f40fdSPierre Jolivet       /* B-part */
700ec4f40fdSPierre Jolivet       for (j=0; j < countB; j++) {
701ec4f40fdSPierre Jolivet         for (k=0; k<bs; k++) {
702ec4f40fdSPierre Jolivet           for (m=0; m<bs; m++) {
703ec4f40fdSPierre Jolivet             if (reuse == MAT_INITIAL_MATRIX) {
7049566063dSJacob Faibussowitsch               PetscCall(PetscMUMPSIntCast(irow + m + shift,&row[jj]));
7059566063dSJacob Faibussowitsch               PetscCall(PetscMUMPSIntCast(garray[bjj[j]]*bs + k + shift,&col[jj]));
706ec4f40fdSPierre Jolivet             }
707ec4f40fdSPierre Jolivet             val[jj++] = v2[j*bs2 + m + k*bs];
708ec4f40fdSPierre Jolivet           }
709ec4f40fdSPierre Jolivet         }
710ec4f40fdSPierre Jolivet       }
711ec4f40fdSPierre Jolivet     } else {
712397b6df1SKris Buschelman       /* A-part */
713397b6df1SKris Buschelman       for (j=0; j<countA; j++) {
714bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
7159566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift,&row[jj]));
7169566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift,&col[jj]));
717397b6df1SKris Buschelman         }
71816ebf90aSShri Abhyankar         val[jj++] = v1[j];
719397b6df1SKris Buschelman       }
72016ebf90aSShri Abhyankar 
72116ebf90aSShri Abhyankar       /* B-part */
72216ebf90aSShri Abhyankar       for (j=0; j < countB; j++) {
723bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
7249566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift,&row[jj]));
7259566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift,&col[jj]));
726397b6df1SKris Buschelman         }
72716ebf90aSShri Abhyankar         val[jj++] = v2[j];
72816ebf90aSShri Abhyankar       }
72916ebf90aSShri Abhyankar     }
730ec4f40fdSPierre Jolivet     irow+=bs;
731ec4f40fdSPierre Jolivet   }
732a6053eceSJunchao Zhang   mumps->nnz = jj;
73316ebf90aSShri Abhyankar   PetscFunctionReturn(0);
73416ebf90aSShri Abhyankar }
73516ebf90aSShri Abhyankar 
736a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
73716ebf90aSShri Abhyankar {
73816ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
739a6053eceSJunchao Zhang   PetscInt64        rstart,nz,i,j,jj,irow,countA,countB;
740a6053eceSJunchao Zhang   PetscMUMPSInt     *row,*col;
74116ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
74216ebf90aSShri Abhyankar   PetscScalar       *val;
743a3d589ffSStefano Zampini   Mat               Ad,Ao;
744a3d589ffSStefano Zampini   Mat_SeqAIJ        *aa;
745a3d589ffSStefano Zampini   Mat_SeqAIJ        *bb;
74616ebf90aSShri Abhyankar 
74716ebf90aSShri Abhyankar   PetscFunctionBegin;
7489566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A,&Ad,&Ao,&garray));
7499566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ad,&av));
7509566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ao,&bv));
751a3d589ffSStefano Zampini 
752a3d589ffSStefano Zampini   aa = (Mat_SeqAIJ*)(Ad)->data;
753a3d589ffSStefano Zampini   bb = (Mat_SeqAIJ*)(Ao)->data;
75438548759SBarry Smith   ai = aa->i;
75538548759SBarry Smith   aj = aa->j;
75638548759SBarry Smith   bi = bb->i;
75738548759SBarry Smith   bj = bb->j;
75816ebf90aSShri Abhyankar 
759a3d589ffSStefano Zampini   rstart = A->rmap->rstart;
7602205254eSKarl Rupp 
761bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
762a6053eceSJunchao Zhang     nz   = (PetscInt64)aa->nz + bb->nz; /* make sure the sum won't overflow PetscInt */
7639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
7649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz,&val));
765a6053eceSJunchao Zhang     mumps->nnz = nz;
766a6053eceSJunchao Zhang     mumps->irn = row;
767a6053eceSJunchao Zhang     mumps->jcn = col;
768a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
76916ebf90aSShri Abhyankar   } else {
770a6053eceSJunchao Zhang     val = mumps->val;
77116ebf90aSShri Abhyankar   }
77216ebf90aSShri Abhyankar 
77316ebf90aSShri Abhyankar   jj = 0; irow = rstart;
77416ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
77516ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
77616ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
77716ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
77816ebf90aSShri Abhyankar     bjj    = bj + bi[i];
77916ebf90aSShri Abhyankar     v1     = av + ai[i];
78016ebf90aSShri Abhyankar     v2     = bv + bi[i];
78116ebf90aSShri Abhyankar 
78216ebf90aSShri Abhyankar     /* A-part */
78316ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
784bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
7859566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift,&row[jj]));
7869566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift,&col[jj]));
78716ebf90aSShri Abhyankar       }
78816ebf90aSShri Abhyankar       val[jj++] = v1[j];
78916ebf90aSShri Abhyankar     }
79016ebf90aSShri Abhyankar 
79116ebf90aSShri Abhyankar     /* B-part */
79216ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
793bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
7949566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift,&row[jj]));
7959566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift,&col[jj]));
79616ebf90aSShri Abhyankar       }
79716ebf90aSShri Abhyankar       val[jj++] = v2[j];
79816ebf90aSShri Abhyankar     }
79916ebf90aSShri Abhyankar     irow++;
80016ebf90aSShri Abhyankar   }
8019566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ad,&av));
8029566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ao,&bv));
80316ebf90aSShri Abhyankar   PetscFunctionReturn(0);
80416ebf90aSShri Abhyankar }
80516ebf90aSShri Abhyankar 
806a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
80767877ebaSShri Abhyankar {
80867877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
80967877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
81067877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
81167877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
812d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
81333d57670SJed Brown   const PetscInt    bs2=mat->bs2;
814a6053eceSJunchao Zhang   PetscInt          bs;
815a6053eceSJunchao Zhang   PetscInt64        nz,i,j,k,n,jj,irow,countA,countB,idx;
816a6053eceSJunchao Zhang   PetscMUMPSInt     *row,*col;
81767877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
81867877ebaSShri Abhyankar   PetscScalar       *val;
81967877ebaSShri Abhyankar 
82067877ebaSShri Abhyankar   PetscFunctionBegin;
8219566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
822bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
82367877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
8249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
8259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz,&val));
826a6053eceSJunchao Zhang     mumps->nnz = nz;
827a6053eceSJunchao Zhang     mumps->irn = row;
828a6053eceSJunchao Zhang     mumps->jcn = col;
829a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
83067877ebaSShri Abhyankar   } else {
831a6053eceSJunchao Zhang     val = mumps->val;
83267877ebaSShri Abhyankar   }
83367877ebaSShri Abhyankar 
834d985c460SShri Abhyankar   jj = 0; irow = rstart;
83567877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
83667877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
83767877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
83867877ebaSShri Abhyankar     ajj    = aj + ai[i];
83967877ebaSShri Abhyankar     bjj    = bj + bi[i];
84067877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
84167877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
84267877ebaSShri Abhyankar 
84367877ebaSShri Abhyankar     idx = 0;
84467877ebaSShri Abhyankar     /* A-part */
84567877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
84667877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
84767877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
848bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
8499566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(irow + n + shift,&row[jj]));
8509566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(rstart + bs*ajj[k] + j + shift,&col[jj]));
85167877ebaSShri Abhyankar           }
85267877ebaSShri Abhyankar           val[jj++] = v1[idx++];
85367877ebaSShri Abhyankar         }
85467877ebaSShri Abhyankar       }
85567877ebaSShri Abhyankar     }
85667877ebaSShri Abhyankar 
85767877ebaSShri Abhyankar     idx = 0;
85867877ebaSShri Abhyankar     /* B-part */
85967877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
86067877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
86167877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
862bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
8639566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(irow + n + shift,&row[jj]));
8649566063dSJacob Faibussowitsch             PetscCall(PetscMUMPSIntCast(bs*garray[bjj[k]] + j + shift,&col[jj]));
86567877ebaSShri Abhyankar           }
866d985c460SShri Abhyankar           val[jj++] = v2[idx++];
86767877ebaSShri Abhyankar         }
86867877ebaSShri Abhyankar       }
86967877ebaSShri Abhyankar     }
870d985c460SShri Abhyankar     irow += bs;
87167877ebaSShri Abhyankar   }
87267877ebaSShri Abhyankar   PetscFunctionReturn(0);
87367877ebaSShri Abhyankar }
87467877ebaSShri Abhyankar 
875a6053eceSJunchao Zhang PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
87616ebf90aSShri Abhyankar {
87716ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
878a6053eceSJunchao Zhang   PetscInt64        rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
879a6053eceSJunchao Zhang   PetscMUMPSInt     *row,*col;
88016ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
88116ebf90aSShri Abhyankar   PetscScalar       *val;
882a3d589ffSStefano Zampini   Mat               Ad,Ao;
883a3d589ffSStefano Zampini   Mat_SeqAIJ        *aa;
884a3d589ffSStefano Zampini   Mat_SeqAIJ        *bb;
88538548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
88638548759SBarry Smith   PetscBool         hermitian;
88738548759SBarry Smith #endif
88816ebf90aSShri Abhyankar 
88916ebf90aSShri Abhyankar   PetscFunctionBegin;
89038548759SBarry Smith #if defined(PETSC_USE_COMPLEX)
8919566063dSJacob Faibussowitsch   PetscCall(MatGetOption(A,MAT_HERMITIAN,&hermitian));
89228b400f6SJacob Faibussowitsch   PetscCheck(!hermitian,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
89338548759SBarry Smith #endif
8949566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A,&Ad,&Ao,&garray));
8959566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ad,&av));
8969566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(Ao,&bv));
897a3d589ffSStefano Zampini 
898a3d589ffSStefano Zampini   aa    = (Mat_SeqAIJ*)(Ad)->data;
899a3d589ffSStefano Zampini   bb    = (Mat_SeqAIJ*)(Ao)->data;
90038548759SBarry Smith   ai    = aa->i;
90138548759SBarry Smith   aj    = aa->j;
90238548759SBarry Smith   adiag = aa->diag;
90338548759SBarry Smith   bi    = bb->i;
90438548759SBarry Smith   bj    = bb->j;
9052205254eSKarl Rupp 
90616ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
90716ebf90aSShri Abhyankar 
908bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
909e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
910e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
91116ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
912e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
91316ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
91416ebf90aSShri Abhyankar       bjj    = bj + bi[i];
915e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
916e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
917e0bace9bSHong Zhang       }
918e0bace9bSHong Zhang     }
91916ebf90aSShri Abhyankar 
920e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
9219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz,&row,nz,&col));
9229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz,&val));
923a6053eceSJunchao Zhang     mumps->nnz = nz;
924a6053eceSJunchao Zhang     mumps->irn = row;
925a6053eceSJunchao Zhang     mumps->jcn = col;
926a6053eceSJunchao Zhang     mumps->val = mumps->val_alloc = val;
92716ebf90aSShri Abhyankar   } else {
928a6053eceSJunchao Zhang     val = mumps->val;
92916ebf90aSShri Abhyankar   }
93016ebf90aSShri Abhyankar 
93116ebf90aSShri Abhyankar   jj = 0; irow = rstart;
93216ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
93316ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
93416ebf90aSShri Abhyankar     v1     = av + adiag[i];
93516ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
93616ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
93716ebf90aSShri Abhyankar     bjj    = bj + bi[i];
93816ebf90aSShri Abhyankar     v2     = bv + bi[i];
93916ebf90aSShri Abhyankar 
94016ebf90aSShri Abhyankar     /* A-part */
94116ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
942bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
9439566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(irow + shift,&row[jj]));
9449566063dSJacob Faibussowitsch         PetscCall(PetscMUMPSIntCast(rstart + ajj[j] + shift,&col[jj]));
94516ebf90aSShri Abhyankar       }
94616ebf90aSShri Abhyankar       val[jj++] = v1[j];
94716ebf90aSShri Abhyankar     }
94816ebf90aSShri Abhyankar 
94916ebf90aSShri Abhyankar     /* B-part */
95016ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
95116ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
952bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
9539566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(irow + shift,&row[jj]));
9549566063dSJacob Faibussowitsch           PetscCall(PetscMUMPSIntCast(garray[bjj[j]] + shift,&col[jj]));
95516ebf90aSShri Abhyankar         }
95616ebf90aSShri Abhyankar         val[jj++] = v2[j];
95716ebf90aSShri Abhyankar       }
958397b6df1SKris Buschelman     }
959397b6df1SKris Buschelman     irow++;
960397b6df1SKris Buschelman   }
9619566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ad,&av));
9629566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(Ao,&bv));
963397b6df1SKris Buschelman   PetscFunctionReturn(0);
964397b6df1SKris Buschelman }
965397b6df1SKris Buschelman 
966dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
967dfbe8321SBarry Smith {
968a6053eceSJunchao Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
969b24902e0SBarry Smith 
970397b6df1SKris Buschelman   PetscFunctionBegin;
9719566063dSJacob Faibussowitsch   PetscCall(PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc));
9729566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mumps->scat_rhs));
9739566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mumps->scat_sol));
9749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mumps->b_seq));
9759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mumps->x_seq));
9769566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.perm_in));
9779566063dSJacob Faibussowitsch   PetscCall(PetscFree2(mumps->irn,mumps->jcn));
9789566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->val_alloc));
9799566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->info));
9809566063dSJacob Faibussowitsch   PetscCall(MatMumpsResetSchur_Private(mumps));
981a5e57a09SHong Zhang   mumps->id.job = JOB_END;
9823ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
983*08401ef6SPierre Jolivet   PetscCheck(mumps->id.INFOG(1) >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in MatDestroy_MUMPS: INFOG(1)=%d",mumps->id.INFOG(1));
9843ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
98567602552SJunchao Zhang   if (mumps->use_petsc_omp_support) {
9869566063dSJacob Faibussowitsch     PetscCall(PetscOmpCtrlDestroy(&mumps->omp_ctrl));
9879566063dSJacob Faibussowitsch     PetscCall(PetscFree2(mumps->rhs_loc,mumps->rhs_recvbuf));
9889566063dSJacob Faibussowitsch     PetscCall(PetscFree3(mumps->rhs_nrow,mumps->rhs_recvcounts,mumps->rhs_disps));
98967602552SJunchao Zhang   }
9903ab56b82SJunchao Zhang #endif
9919566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->ia_alloc));
9929566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->ja_alloc));
9939566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->recvcount));
9949566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->reqs));
9959566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->irhs_loc));
9969566063dSJacob Faibussowitsch   if (mumps->mumps_comm != MPI_COMM_NULL) PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)A),&mumps->mumps_comm));
9979566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
998bf0cc555SLisandro Dalcin 
99997969023SHong Zhang   /* clear composed functions */
10009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL));
10019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL));
10029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL));
10039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL));
10049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL));
10059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL));
10069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL));
10079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL));
10089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL));
10099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL));
10109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL));
10119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL));
10129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverseTranspose_C",NULL));
1013397b6df1SKris Buschelman   PetscFunctionReturn(0);
1014397b6df1SKris Buschelman }
1015397b6df1SKris Buschelman 
101667602552SJunchao 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. */
101767602552SJunchao Zhang static PetscErrorCode MatMumpsSetUpDistRHSInfo(Mat A,PetscInt nrhs,const PetscScalar *array)
101867602552SJunchao Zhang {
101967602552SJunchao Zhang   Mat_MUMPS          *mumps=(Mat_MUMPS*)A->data;
102067602552SJunchao Zhang   const PetscMPIInt  ompsize=mumps->omp_comm_size;
102167602552SJunchao Zhang   PetscInt           i,m,M,rstart;
102267602552SJunchao Zhang 
102367602552SJunchao Zhang   PetscFunctionBegin;
10249566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,NULL));
10259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,NULL));
1026*08401ef6SPierre Jolivet   PetscCheck(M <= PETSC_MUMPS_INT_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PetscInt too long for PetscMUMPSInt");
102767602552SJunchao Zhang   if (ompsize == 1) {
102867602552SJunchao Zhang     if (!mumps->irhs_loc) {
102967602552SJunchao Zhang       mumps->nloc_rhs = m;
10309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m,&mumps->irhs_loc));
10319566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(A,&rstart,NULL));
103267602552SJunchao Zhang       for (i=0; i<m; i++) mumps->irhs_loc[i] = rstart+i+1; /* use 1-based indices */
103367602552SJunchao Zhang     }
103467602552SJunchao Zhang     mumps->id.rhs_loc = (MumpsScalar*)array;
103567602552SJunchao Zhang   } else {
103667602552SJunchao Zhang   #if defined(PETSC_HAVE_OPENMP_SUPPORT)
103767602552SJunchao Zhang     const PetscInt  *ranges;
103867602552SJunchao Zhang     PetscMPIInt     j,k,sendcount,*petsc_ranks,*omp_ranks;
103967602552SJunchao Zhang     MPI_Group       petsc_group,omp_group;
104067602552SJunchao Zhang     PetscScalar     *recvbuf=NULL;
104167602552SJunchao Zhang 
104267602552SJunchao Zhang     if (mumps->is_omp_master) {
104367602552SJunchao Zhang       /* Lazily initialize the omp stuff for distributed rhs */
104467602552SJunchao Zhang       if (!mumps->irhs_loc) {
10459566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(ompsize,&omp_ranks,ompsize,&petsc_ranks));
10469566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(ompsize,&mumps->rhs_nrow,ompsize,&mumps->rhs_recvcounts,ompsize,&mumps->rhs_disps));
10479566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_group(mumps->petsc_comm,&petsc_group));
10489566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_group(mumps->omp_comm,&omp_group));
104967602552SJunchao Zhang         for (j=0; j<ompsize; j++) omp_ranks[j] = j;
10509566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_translate_ranks(omp_group,ompsize,omp_ranks,petsc_group,petsc_ranks));
105167602552SJunchao Zhang 
105267602552SJunchao Zhang         /* Populate mumps->irhs_loc[], rhs_nrow[] */
105367602552SJunchao Zhang         mumps->nloc_rhs = 0;
10549566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRanges(A,&ranges));
105567602552SJunchao Zhang         for (j=0; j<ompsize; j++) {
105667602552SJunchao Zhang           mumps->rhs_nrow[j] = ranges[petsc_ranks[j]+1] - ranges[petsc_ranks[j]];
105767602552SJunchao Zhang           mumps->nloc_rhs   += mumps->rhs_nrow[j];
105867602552SJunchao Zhang         }
10599566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(mumps->nloc_rhs,&mumps->irhs_loc));
106067602552SJunchao Zhang         for (j=k=0; j<ompsize; j++) {
106167602552SJunchao 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 */
106267602552SJunchao Zhang         }
106367602552SJunchao Zhang 
10649566063dSJacob Faibussowitsch         PetscCall(PetscFree2(omp_ranks,petsc_ranks));
10659566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_free(&petsc_group));
10669566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Group_free(&omp_group));
106767602552SJunchao Zhang       }
106867602552SJunchao Zhang 
106967602552SJunchao Zhang       /* Realloc buffers when current nrhs is bigger than what we have met */
107067602552SJunchao Zhang       if (nrhs > mumps->max_nrhs) {
10719566063dSJacob Faibussowitsch         PetscCall(PetscFree2(mumps->rhs_loc,mumps->rhs_recvbuf));
10729566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(mumps->nloc_rhs*nrhs,&mumps->rhs_loc,mumps->nloc_rhs*nrhs,&mumps->rhs_recvbuf));
107367602552SJunchao Zhang         mumps->max_nrhs = nrhs;
107467602552SJunchao Zhang       }
107567602552SJunchao Zhang 
107667602552SJunchao Zhang       /* Setup recvcounts[], disps[], recvbuf on omp rank 0 for the upcoming MPI_Gatherv */
10779566063dSJacob Faibussowitsch       for (j=0; j<ompsize; j++) PetscCall(PetscMPIIntCast(mumps->rhs_nrow[j]*nrhs,&mumps->rhs_recvcounts[j]));
107867602552SJunchao Zhang       mumps->rhs_disps[0] = 0;
107967602552SJunchao Zhang       for (j=1; j<ompsize; j++) {
108067602552SJunchao Zhang         mumps->rhs_disps[j] = mumps->rhs_disps[j-1] + mumps->rhs_recvcounts[j-1];
1081*08401ef6SPierre Jolivet         PetscCheck(mumps->rhs_disps[j] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PetscMPIInt overflow!");
108267602552SJunchao Zhang       }
108367602552SJunchao Zhang       recvbuf = (nrhs == 1) ? mumps->rhs_loc : mumps->rhs_recvbuf; /* Directly use rhs_loc[] as recvbuf. Single rhs is common in Ax=b */
108467602552SJunchao Zhang     }
108567602552SJunchao Zhang 
10869566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(m*nrhs,&sendcount));
10879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(array,sendcount,MPIU_SCALAR,recvbuf,mumps->rhs_recvcounts,mumps->rhs_disps,MPIU_SCALAR,0,mumps->omp_comm));
108867602552SJunchao Zhang 
108967602552SJunchao Zhang     if (mumps->is_omp_master) {
109067602552SJunchao Zhang       if (nrhs > 1) { /* Copy & re-arrange data from rhs_recvbuf[] to mumps->rhs_loc[] only when there are multiple rhs */
109167602552SJunchao Zhang         PetscScalar *dst,*dstbase = mumps->rhs_loc;
109267602552SJunchao Zhang         for (j=0; j<ompsize; j++) {
109367602552SJunchao Zhang           const PetscScalar *src = mumps->rhs_recvbuf + mumps->rhs_disps[j];
109467602552SJunchao Zhang           dst = dstbase;
109567602552SJunchao Zhang           for (i=0; i<nrhs; i++) {
10969566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(dst,src,mumps->rhs_nrow[j]));
109767602552SJunchao Zhang             src += mumps->rhs_nrow[j];
109867602552SJunchao Zhang             dst += mumps->nloc_rhs;
109967602552SJunchao Zhang           }
110067602552SJunchao Zhang           dstbase += mumps->rhs_nrow[j];
110167602552SJunchao Zhang         }
110267602552SJunchao Zhang       }
110367602552SJunchao Zhang       mumps->id.rhs_loc = (MumpsScalar*)mumps->rhs_loc;
110467602552SJunchao Zhang     }
110567602552SJunchao Zhang   #endif /* PETSC_HAVE_OPENMP_SUPPORT */
110667602552SJunchao Zhang   }
110767602552SJunchao Zhang   mumps->id.nrhs     = nrhs;
110867602552SJunchao Zhang   mumps->id.nloc_rhs = mumps->nloc_rhs;
110967602552SJunchao Zhang   mumps->id.lrhs_loc = mumps->nloc_rhs;
111067602552SJunchao Zhang   mumps->id.irhs_loc = mumps->irhs_loc;
111167602552SJunchao Zhang   PetscFunctionReturn(0);
111267602552SJunchao Zhang }
111367602552SJunchao Zhang 
1114b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
1115b24902e0SBarry Smith {
1116e69c285eSBarry Smith   Mat_MUMPS          *mumps=(Mat_MUMPS*)A->data;
111725aac85cSJunchao Zhang   const PetscScalar  *rarray = NULL;
1118d54de34fSKris Buschelman   PetscScalar        *array;
1119329ec9b3SHong Zhang   IS                 is_iden,is_petsc;
1120329ec9b3SHong Zhang   PetscInt           i;
1121cc86f929SStefano Zampini   PetscBool          second_solve = PETSC_FALSE;
1122883f2eb9SBarry Smith   static PetscBool   cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
1123397b6df1SKris Buschelman 
1124397b6df1SKris Buschelman   PetscFunctionBegin;
11259566063dSJacob Faibussowitsch   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 Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",&cite1));
11269566063dSJacob Faibussowitsch   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 Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",&cite2));
11272aca8efcSHong Zhang 
1128603e8f96SBarry Smith   if (A->factorerrortype) {
11299566063dSJacob 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)));
11309566063dSJacob Faibussowitsch     PetscCall(VecSetInf(x));
11312aca8efcSHong Zhang     PetscFunctionReturn(0);
11322aca8efcSHong Zhang   }
11332aca8efcSHong Zhang 
1134a5e57a09SHong Zhang   mumps->id.nrhs = 1;
11352d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
113625aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
113767602552SJunchao Zhang       mumps->id.ICNTL(20) = 10; /* dense distributed RHS */
11389566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(b,&rarray));
11399566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetUpDistRHSInfo(A,1,rarray));
114025aac85cSJunchao Zhang     } else {
114141ffd417SStefano Zampini       mumps->id.ICNTL(20) = 0; /* dense centralized RHS; Scatter b into a sequential rhs vector*/
11429566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mumps->scat_rhs,b,mumps->b_seq,INSERT_VALUES,SCATTER_FORWARD));
11439566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mumps->scat_rhs,b,mumps->b_seq,INSERT_VALUES,SCATTER_FORWARD));
114467602552SJunchao Zhang       if (!mumps->myid) {
11459566063dSJacob Faibussowitsch         PetscCall(VecGetArray(mumps->b_seq,&array));
114667602552SJunchao Zhang         mumps->id.rhs = (MumpsScalar*)array;
114767602552SJunchao Zhang       }
114825aac85cSJunchao Zhang     }
11493ab56b82SJunchao Zhang   } else {  /* petsc_size == 1 */
115067602552SJunchao Zhang     mumps->id.ICNTL(20) = 0; /* dense centralized RHS */
11519566063dSJacob Faibussowitsch     PetscCall(VecCopy(b,x));
11529566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x,&array));
1153940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
1154397b6df1SKris Buschelman   }
1155397b6df1SKris Buschelman 
1156cc86f929SStefano Zampini   /*
1157cc86f929SStefano Zampini      handle condensation step of Schur complement (if any)
1158cc86f929SStefano Zampini      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
1159cc86f929SStefano Zampini      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
1160cc86f929SStefano Zampini      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
1161cc86f929SStefano Zampini      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
1162cc86f929SStefano Zampini   */
1163583f777eSStefano Zampini   if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
1164*08401ef6SPierre Jolivet     PetscCheck(mumps->petsc_size <= 1,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc");
1165cc86f929SStefano Zampini     second_solve = PETSC_TRUE;
11669566063dSJacob Faibussowitsch     PetscCall(MatMumpsHandleSchur_Private(A,PETSC_FALSE));
1167cc86f929SStefano Zampini   }
1168397b6df1SKris Buschelman   /* solve phase */
1169329ec9b3SHong Zhang   /*-------------*/
1170a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
11713ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
1172*08401ef6SPierre Jolivet   PetscCheck(mumps->id.INFOG(1) >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d",mumps->id.INFOG(1));
1173397b6df1SKris Buschelman 
1174b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
1175cc86f929SStefano Zampini   if (second_solve) {
11769566063dSJacob Faibussowitsch     PetscCall(MatMumpsHandleSchur_Private(A,PETSC_TRUE));
1177cc86f929SStefano Zampini   }
1178b5fa320bSStefano Zampini 
11792d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) { /* convert mumps distributed solution to petsc mpi x */
1180a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
1181a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
11829566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&mumps->scat_sol));
1183397b6df1SKris Buschelman     }
1184a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
1185a6053eceSJunchao Zhang       PetscInt *isol2_loc=NULL;
11869566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden)); /* from */
11879566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mumps->id.lsol_loc,&isol2_loc));
1188a6053eceSJunchao 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 */
11899566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,isol2_loc,PETSC_OWN_POINTER,&is_petsc));  /* to */
11909566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol));
11919566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_iden));
11929566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_petsc));
1193a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
1194397b6df1SKris Buschelman     }
1195a5e57a09SHong Zhang 
11969566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD));
11979566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD));
1198329ec9b3SHong Zhang   }
1199353d7d71SJunchao Zhang 
120067602552SJunchao Zhang   if (mumps->petsc_size > 1) {
120125aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
12029566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(b,&rarray));
120325aac85cSJunchao Zhang     } else if (!mumps->myid) {
12049566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mumps->b_seq,&array));
120525aac85cSJunchao Zhang     }
12069566063dSJacob Faibussowitsch   } else PetscCall(VecRestoreArray(x,&array));
1207353d7d71SJunchao Zhang 
12089566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0*mumps->id.RINFO(3)));
1209397b6df1SKris Buschelman   PetscFunctionReturn(0);
1210397b6df1SKris Buschelman }
1211397b6df1SKris Buschelman 
121251d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
121351d5961aSHong Zhang {
1214e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
121551d5961aSHong Zhang 
121651d5961aSHong Zhang   PetscFunctionBegin;
1217a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
12189566063dSJacob Faibussowitsch   PetscCall(MatSolve_MUMPS(A,b,x));
1219a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
122051d5961aSHong Zhang   PetscFunctionReturn(0);
122151d5961aSHong Zhang }
122251d5961aSHong Zhang 
1223e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
1224e0b74bf9SHong Zhang {
1225b8491c3eSStefano Zampini   Mat               Bt = NULL;
1226a6053eceSJunchao Zhang   PetscBool         denseX,denseB,flg,flgT;
1227e69c285eSBarry Smith   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
1228334c5f61SHong Zhang   PetscInt          i,nrhs,M;
12291683a169SBarry Smith   PetscScalar       *array;
12301683a169SBarry Smith   const PetscScalar *rbray;
1231a6053eceSJunchao Zhang   PetscInt          lsol_loc,nlsol_loc,*idxx,iidx = 0;
1232a6053eceSJunchao Zhang   PetscMUMPSInt     *isol_loc,*isol_loc_save;
12331683a169SBarry Smith   PetscScalar       *bray,*sol_loc,*sol_loc_save;
1234be818407SHong Zhang   IS                is_to,is_from;
1235beae5ec0SHong Zhang   PetscInt          k,proc,j,m,myrstart;
1236be818407SHong Zhang   const PetscInt    *rstart;
123767602552SJunchao Zhang   Vec               v_mpi,msol_loc;
123867602552SJunchao Zhang   VecScatter        scat_sol;
123967602552SJunchao Zhang   Vec               b_seq;
124067602552SJunchao Zhang   VecScatter        scat_rhs;
1241be818407SHong Zhang   PetscScalar       *aa;
1242be818407SHong Zhang   PetscInt          spnr,*ia,*ja;
1243d56c302dSHong Zhang   Mat_MPIAIJ        *b = NULL;
1244bda8bf91SBarry Smith 
1245e0b74bf9SHong Zhang   PetscFunctionBegin;
12469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X,&denseX,MATSEQDENSE,MATMPIDENSE,NULL));
124728b400f6SJacob Faibussowitsch   PetscCheck(denseX,PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
1248be818407SHong Zhang 
12499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&denseB,MATSEQDENSE,MATMPIDENSE,NULL));
1250a6053eceSJunchao Zhang   if (denseB) {
1251*08401ef6SPierre Jolivet     PetscCheck(B->rmap->n == X->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
1252be818407SHong Zhang     mumps->id.ICNTL(20)= 0; /* dense RHS */
12530e6b8875SHong Zhang   } else { /* sparse B */
1254*08401ef6SPierre Jolivet     PetscCheck(X != B,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_IDN,"X and B must be different matrices");
12559566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT));
12560e6b8875SHong Zhang     if (flgT) { /* input B is transpose of actural RHS matrix,
12570e6b8875SHong Zhang                  because mumps requires sparse compressed COLUMN storage! See MatMatTransposeSolve_MUMPS() */
12589566063dSJacob Faibussowitsch       PetscCall(MatTransposeGetMat(B,&Bt));
12590f52d626SJunchao Zhang     } else SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATTRANSPOSEMAT matrix");
1260be818407SHong Zhang     mumps->id.ICNTL(20)= 1; /* sparse RHS */
1261b8491c3eSStefano Zampini   }
126287b22cf4SHong Zhang 
12639566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B,&M,&nrhs));
12649481e6e9SHong Zhang   mumps->id.nrhs = nrhs;
12659481e6e9SHong Zhang   mumps->id.lrhs = M;
12662b691707SHong Zhang   mumps->id.rhs  = NULL;
12679481e6e9SHong Zhang 
12682d4298aeSJunchao Zhang   if (mumps->petsc_size == 1) {
1269b8491c3eSStefano Zampini     PetscScalar *aa;
1270b8491c3eSStefano Zampini     PetscInt    spnr,*ia,*ja;
1271e94cce23SStefano Zampini     PetscBool   second_solve = PETSC_FALSE;
1272b8491c3eSStefano Zampini 
12739566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(X,&array));
1274b8491c3eSStefano Zampini     mumps->id.rhs = (MumpsScalar*)array;
12752b691707SHong Zhang 
1276a6053eceSJunchao Zhang     if (denseB) {
12772b691707SHong Zhang       /* copy B to X */
12789566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(B,&rbray));
12799566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array,rbray,M*nrhs));
12809566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(B,&rbray));
12812b691707SHong Zhang     } else { /* sparse B */
12829566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(Bt,&aa));
12839566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg));
128428b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
12859566063dSJacob Faibussowitsch       PetscCall(PetscMUMPSIntCSRCast(mumps,spnr,ia,ja,&mumps->id.irhs_ptr,&mumps->id.irhs_sparse,&mumps->id.nz_rhs));
1286b8491c3eSStefano Zampini       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1287b8491c3eSStefano Zampini     }
1288e94cce23SStefano Zampini     /* handle condensation step of Schur complement (if any) */
1289583f777eSStefano Zampini     if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
1290e94cce23SStefano Zampini       second_solve = PETSC_TRUE;
12919566063dSJacob Faibussowitsch       PetscCall(MatMumpsHandleSchur_Private(A,PETSC_FALSE));
1292e94cce23SStefano Zampini     }
12932cd7d884SHong Zhang     /* solve phase */
12942cd7d884SHong Zhang     /*-------------*/
12952cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
12963ab56b82SJunchao Zhang     PetscMUMPS_c(mumps);
1297*08401ef6SPierre Jolivet     PetscCheck(mumps->id.INFOG(1) >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d",mumps->id.INFOG(1));
1298b5fa320bSStefano Zampini 
1299b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
1300e94cce23SStefano Zampini     if (second_solve) {
13019566063dSJacob Faibussowitsch       PetscCall(MatMumpsHandleSchur_Private(A,PETSC_TRUE));
1302e94cce23SStefano Zampini     }
1303a6053eceSJunchao Zhang     if (!denseB) { /* sparse B */
13049566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(Bt,&aa));
13059566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg));
130628b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
1307b8491c3eSStefano Zampini     }
13089566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X,&array));
1309be818407SHong Zhang     PetscFunctionReturn(0);
1310be818407SHong Zhang   }
1311801fbe65SHong Zhang 
1312be818407SHong Zhang   /*--------- parallel case: MUMPS requires rhs B to be centralized on the host! --------*/
13132c71b3e2SJacob Faibussowitsch   PetscCheckFalse(mumps->petsc_size > 1 && mumps->id.ICNTL(19),PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc");
1314241dbb5eSStefano Zampini 
1315beae5ec0SHong Zhang   /* create msol_loc to hold mumps local solution */
13161683a169SBarry Smith   isol_loc_save = mumps->id.isol_loc; /* save it for MatSolve() */
13171683a169SBarry Smith   sol_loc_save  = (PetscScalar*)mumps->id.sol_loc;
1318801fbe65SHong Zhang 
1319a1dfcbd9SJunchao Zhang   lsol_loc  = mumps->id.lsol_loc;
132071aed81dSHong Zhang   nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
13219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nlsol_loc,&sol_loc,lsol_loc,&isol_loc));
1322940cd9d6SSatish Balay   mumps->id.sol_loc  = (MumpsScalar*)sol_loc;
1323801fbe65SHong Zhang   mumps->id.isol_loc = isol_loc;
1324801fbe65SHong Zhang 
13259566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&msol_loc));
13262cd7d884SHong Zhang 
132767602552SJunchao Zhang   if (denseB) {
132825aac85cSJunchao Zhang     if (mumps->ICNTL20 == 10) {
132967602552SJunchao Zhang       mumps->id.ICNTL(20) = 10; /* dense distributed RHS */
13309566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(B,&rbray));
13319566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetUpDistRHSInfo(A,nrhs,rbray));
13329566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(B,&rbray));
13339566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B,&m,NULL));
13349566063dSJacob Faibussowitsch       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,NULL,&v_mpi));
133525aac85cSJunchao Zhang     } else {
133625aac85cSJunchao Zhang       mumps->id.ICNTL(20) = 0; /* dense centralized RHS */
133780577c12SJunchao Zhang       /* TODO: Because of non-contiguous indices, the created vecscatter scat_rhs is not done in MPI_Gather, resulting in
133880577c12SJunchao Zhang         very inefficient communication. An optimization is to use VecScatterCreateToZero to gather B to rank 0. Then on rank
133980577c12SJunchao Zhang         0, re-arrange B into desired order, which is a local operation.
134080577c12SJunchao Zhang       */
134180577c12SJunchao Zhang 
134267602552SJunchao Zhang       /* scatter v_mpi to b_seq because MUMPS before 5.3.0 only supports centralized rhs */
1343be818407SHong Zhang       /* wrap dense rhs matrix B into a vector v_mpi */
13449566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B,&m,NULL));
13459566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArray(B,&bray));
13469566063dSJacob Faibussowitsch       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi));
13479566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArray(B,&bray));
13482b691707SHong Zhang 
1349be818407SHong Zhang       /* scatter v_mpi to b_seq in proc[0]. MUMPS requires rhs to be centralized on the host! */
1350801fbe65SHong Zhang       if (!mumps->myid) {
1351beae5ec0SHong Zhang         PetscInt *idx;
1352beae5ec0SHong Zhang         /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B */
13539566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nrhs*M,&idx));
13549566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRanges(B,&rstart));
1355be818407SHong Zhang         k = 0;
13562d4298aeSJunchao Zhang         for (proc=0; proc<mumps->petsc_size; proc++) {
1357be818407SHong Zhang           for (j=0; j<nrhs; j++) {
1358beae5ec0SHong Zhang             for (i=rstart[proc]; i<rstart[proc+1]; i++) idx[k++] = j*M + i;
1359be818407SHong Zhang           }
1360be818407SHong Zhang         }
1361be818407SHong Zhang 
13629566063dSJacob Faibussowitsch         PetscCall(VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq));
13639566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_OWN_POINTER,&is_to));
13649566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from));
1365801fbe65SHong Zhang       } else {
13669566063dSJacob Faibussowitsch         PetscCall(VecCreateSeq(PETSC_COMM_SELF,0,&b_seq));
13679566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to));
13689566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from));
1369801fbe65SHong Zhang       }
13709566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs));
13719566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD));
13729566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_to));
13739566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_from));
13749566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD));
1375801fbe65SHong Zhang 
1376801fbe65SHong Zhang       if (!mumps->myid) { /* define rhs on the host */
13779566063dSJacob Faibussowitsch         PetscCall(VecGetArray(b_seq,&bray));
1378940cd9d6SSatish Balay         mumps->id.rhs = (MumpsScalar*)bray;
13799566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(b_seq,&bray));
1380801fbe65SHong Zhang       }
138125aac85cSJunchao Zhang     }
13822b691707SHong Zhang   } else { /* sparse B */
13832b691707SHong Zhang     b = (Mat_MPIAIJ*)Bt->data;
13842b691707SHong Zhang 
1385be818407SHong Zhang     /* wrap dense X into a vector v_mpi */
13869566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(X,&m,NULL));
13879566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(X,&bray));
13889566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi));
13899566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X,&bray));
13902b691707SHong Zhang 
13912b691707SHong Zhang     if (!mumps->myid) {
13929566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(b->A,&aa));
13939566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg));
139428b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
13959566063dSJacob Faibussowitsch       PetscCall(PetscMUMPSIntCSRCast(mumps,spnr,ia,ja,&mumps->id.irhs_ptr,&mumps->id.irhs_sparse,&mumps->id.nz_rhs));
13962b691707SHong Zhang       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
13972b691707SHong Zhang     } else {
13982b691707SHong Zhang       mumps->id.irhs_ptr    = NULL;
13992b691707SHong Zhang       mumps->id.irhs_sparse = NULL;
14002b691707SHong Zhang       mumps->id.nz_rhs      = 0;
14012b691707SHong Zhang       mumps->id.rhs_sparse  = NULL;
14022b691707SHong Zhang     }
14032b691707SHong Zhang   }
14042b691707SHong Zhang 
1405801fbe65SHong Zhang   /* solve phase */
1406801fbe65SHong Zhang   /*-------------*/
1407801fbe65SHong Zhang   mumps->id.job = JOB_SOLVE;
14083ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
1409*08401ef6SPierre Jolivet   PetscCheck(mumps->id.INFOG(1) >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d",mumps->id.INFOG(1));
1410801fbe65SHong Zhang 
1411334c5f61SHong Zhang   /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
14129566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(X,&array));
14139566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(v_mpi,array));
1414801fbe65SHong Zhang 
1415334c5f61SHong Zhang   /* create scatter scat_sol */
14169566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(X,&rstart));
1417beae5ec0SHong Zhang   /* iidx: index for scatter mumps solution to petsc X */
1418beae5ec0SHong Zhang 
14199566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from));
14209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nlsol_loc,&idxx));
1421beae5ec0SHong Zhang   for (i=0; i<lsol_loc; i++) {
1422beae5ec0SHong 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 */
1423beae5ec0SHong Zhang 
14242d4298aeSJunchao Zhang     for (proc=0; proc<mumps->petsc_size; proc++) {
1425beae5ec0SHong Zhang       if (isol_loc[i] >= rstart[proc] && isol_loc[i] < rstart[proc+1]) {
1426beae5ec0SHong Zhang         myrstart = rstart[proc];
1427beae5ec0SHong Zhang         k        = isol_loc[i] - myrstart;        /* local index on 1st column of petsc vector X */
1428beae5ec0SHong Zhang         iidx     = k + myrstart*nrhs;             /* maps mumps isol_loc[i] to petsc index in X */
1429beae5ec0SHong Zhang         m        = rstart[proc+1] - rstart[proc]; /* rows of X for this proc */
1430beae5ec0SHong Zhang         break;
1431be818407SHong Zhang       }
1432be818407SHong Zhang     }
1433be818407SHong Zhang 
1434beae5ec0SHong Zhang     for (j=0; j<nrhs; j++) idxx[i+j*lsol_loc] = iidx + j*m;
1435801fbe65SHong Zhang   }
14369566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to));
14379566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(msol_loc,is_from,v_mpi,is_to,&scat_sol));
14389566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat_sol,msol_loc,v_mpi,INSERT_VALUES,SCATTER_FORWARD));
14399566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_from));
14409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_to));
14419566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat_sol,msol_loc,v_mpi,INSERT_VALUES,SCATTER_FORWARD));
14429566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(X,&array));
144371aed81dSHong Zhang 
144471aed81dSHong Zhang   /* free spaces */
14451683a169SBarry Smith   mumps->id.sol_loc  = (MumpsScalar*)sol_loc_save;
144671aed81dSHong Zhang   mumps->id.isol_loc = isol_loc_save;
144771aed81dSHong Zhang 
14489566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sol_loc,isol_loc));
14499566063dSJacob Faibussowitsch   PetscCall(PetscFree(idxx));
14509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&msol_loc));
14519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v_mpi));
1452a6053eceSJunchao Zhang   if (!denseB) {
14532b691707SHong Zhang     if (!mumps->myid) {
1454d56c302dSHong Zhang       b = (Mat_MPIAIJ*)Bt->data;
14559566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(b->A,&aa));
14569566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg));
145728b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
14582b691707SHong Zhang     }
14592b691707SHong Zhang   } else {
146025aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) {
14619566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b_seq));
14629566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&scat_rhs));
146325aac85cSJunchao Zhang     }
14642b691707SHong Zhang   }
14659566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scat_sol));
14669566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0*nrhs*mumps->id.RINFO(3)));
1467e0b74bf9SHong Zhang   PetscFunctionReturn(0);
1468e0b74bf9SHong Zhang }
1469e0b74bf9SHong Zhang 
1470eb3ef3b2SHong Zhang PetscErrorCode MatMatTransposeSolve_MUMPS(Mat A,Mat Bt,Mat X)
1471eb3ef3b2SHong Zhang {
1472eb3ef3b2SHong Zhang   PetscBool      flg;
1473eb3ef3b2SHong Zhang   Mat            B;
1474eb3ef3b2SHong Zhang 
1475eb3ef3b2SHong Zhang   PetscFunctionBegin;
14769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)Bt,&flg,MATSEQAIJ,MATMPIAIJ,NULL));
147728b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Matrix Bt must be MATAIJ matrix");
1478eb3ef3b2SHong Zhang 
1479eb3ef3b2SHong Zhang   /* Create B=Bt^T that uses Bt's data structure */
14809566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(Bt,&B));
1481eb3ef3b2SHong Zhang 
14829566063dSJacob Faibussowitsch   PetscCall(MatMatSolve_MUMPS(A,B,X));
14839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1484eb3ef3b2SHong Zhang   PetscFunctionReturn(0);
1485eb3ef3b2SHong Zhang }
1486eb3ef3b2SHong Zhang 
1487ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1488a58c3f20SHong Zhang /*
1489a58c3f20SHong Zhang   input:
1490a58c3f20SHong Zhang    F:        numeric factor
1491a58c3f20SHong Zhang   output:
1492a58c3f20SHong Zhang    nneg:     total number of negative pivots
149319d49a3bSHong Zhang    nzero:    total number of zero pivots
149419d49a3bSHong Zhang    npos:     (global dimension of F) - nneg - nzero
1495a58c3f20SHong Zhang */
1496a6053eceSJunchao Zhang PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
1497a58c3f20SHong Zhang {
1498e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1499c1490034SHong Zhang   PetscMPIInt    size;
1500a58c3f20SHong Zhang 
1501a58c3f20SHong Zhang   PetscFunctionBegin;
15029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)F),&size));
1503bcb30aebSHong 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 */
15042c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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));
1505ed85ac9fSHong Zhang 
1506710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1507ed85ac9fSHong Zhang   if (nzero || npos) {
1508*08401ef6SPierre 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");
1509710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1510710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1511a58c3f20SHong Zhang   }
1512a58c3f20SHong Zhang   PetscFunctionReturn(0);
1513a58c3f20SHong Zhang }
151419d49a3bSHong Zhang #endif
1515a58c3f20SHong Zhang 
15163ab56b82SJunchao Zhang PetscErrorCode MatMumpsGatherNonzerosOnMaster(MatReuse reuse,Mat_MUMPS *mumps)
15173ab56b82SJunchao Zhang {
1518a6053eceSJunchao Zhang   PetscInt       i,nreqs;
1519a6053eceSJunchao Zhang   PetscMUMPSInt  *irn,*jcn;
1520a6053eceSJunchao Zhang   PetscMPIInt    count;
1521a6053eceSJunchao Zhang   PetscInt64     totnnz,remain;
1522a6053eceSJunchao Zhang   const PetscInt osize=mumps->omp_comm_size;
1523a6053eceSJunchao Zhang   PetscScalar    *val;
15243ab56b82SJunchao Zhang 
15253ab56b82SJunchao Zhang   PetscFunctionBegin;
1526a6053eceSJunchao Zhang   if (osize > 1) {
15273ab56b82SJunchao Zhang     if (reuse == MAT_INITIAL_MATRIX) {
15283ab56b82SJunchao Zhang       /* master first gathers counts of nonzeros to receive */
15299566063dSJacob Faibussowitsch       if (mumps->is_omp_master) PetscCall(PetscMalloc1(osize,&mumps->recvcount));
15309566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gather(&mumps->nnz,1,MPIU_INT64,mumps->recvcount,1,MPIU_INT64,0/*master*/,mumps->omp_comm));
15313ab56b82SJunchao Zhang 
1532a6053eceSJunchao Zhang       /* Then each computes number of send/recvs */
15333ab56b82SJunchao Zhang       if (mumps->is_omp_master) {
1534a6053eceSJunchao Zhang         /* Start from 1 since self communication is not done in MPI */
1535a6053eceSJunchao Zhang         nreqs = 0;
1536a6053eceSJunchao Zhang         for (i=1; i<osize; i++) nreqs += (mumps->recvcount[i]+PETSC_MPI_INT_MAX-1)/PETSC_MPI_INT_MAX;
1537a6053eceSJunchao Zhang       } else {
1538a6053eceSJunchao Zhang         nreqs = (mumps->nnz+PETSC_MPI_INT_MAX-1)/PETSC_MPI_INT_MAX;
15393ab56b82SJunchao Zhang       }
15409566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nreqs*3,&mumps->reqs)); /* Triple the requests since we send irn, jcn and val seperately */
15413ab56b82SJunchao Zhang 
1542a6053eceSJunchao Zhang       /* The following code is doing a very simple thing: omp_master rank gathers irn/jcn/val from others.
1543a6053eceSJunchao Zhang          MPI_Gatherv would be enough if it supports big counts > 2^31-1. Since it does not, and mumps->nnz
1544a6053eceSJunchao Zhang          might be a prime number > 2^31-1, we have to slice the message. Note omp_comm_size
1545a6053eceSJunchao Zhang          is very small, the current approach should have no extra overhead compared to MPI_Gatherv.
1546a6053eceSJunchao Zhang        */
1547a6053eceSJunchao Zhang       nreqs = 0; /* counter for actual send/recvs */
15483ab56b82SJunchao Zhang       if (mumps->is_omp_master) {
1549a6053eceSJunchao Zhang         for (i=0,totnnz=0; i<osize; i++) totnnz += mumps->recvcount[i]; /* totnnz = sum of nnz over omp_comm */
15509566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(totnnz,&irn,totnnz,&jcn));
15519566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(totnnz,&val));
1552a6053eceSJunchao Zhang 
1553a6053eceSJunchao Zhang         /* Self communication */
15549566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(irn,mumps->irn,mumps->nnz));
15559566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jcn,mumps->jcn,mumps->nnz));
15569566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(val,mumps->val,mumps->nnz));
1557a6053eceSJunchao Zhang 
1558a6053eceSJunchao Zhang         /* Replace mumps->irn/jcn etc on master with the newly allocated bigger arrays */
15599566063dSJacob Faibussowitsch         PetscCall(PetscFree2(mumps->irn,mumps->jcn));
15609566063dSJacob Faibussowitsch         PetscCall(PetscFree(mumps->val_alloc));
1561a6053eceSJunchao Zhang         mumps->nnz = totnnz;
15623ab56b82SJunchao Zhang         mumps->irn = irn;
15633ab56b82SJunchao Zhang         mumps->jcn = jcn;
1564a6053eceSJunchao Zhang         mumps->val = mumps->val_alloc = val;
1565a6053eceSJunchao Zhang 
1566a6053eceSJunchao Zhang         irn += mumps->recvcount[0]; /* recvcount[0] is old mumps->nnz on omp rank 0 */
1567a6053eceSJunchao Zhang         jcn += mumps->recvcount[0];
1568a6053eceSJunchao Zhang         val += mumps->recvcount[0];
1569a6053eceSJunchao Zhang 
1570a6053eceSJunchao Zhang         /* Remote communication */
1571a6053eceSJunchao Zhang         for (i=1; i<osize; i++) {
1572a6053eceSJunchao Zhang           count  = PetscMin(mumps->recvcount[i],PETSC_MPI_INT_MAX);
1573a6053eceSJunchao Zhang           remain = mumps->recvcount[i] - count;
1574a6053eceSJunchao Zhang           while (count>0) {
15759566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(irn,count,MPIU_MUMPSINT,i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
15769566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(jcn,count,MPIU_MUMPSINT,i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
15779566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(val,count,MPIU_SCALAR,  i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
1578a6053eceSJunchao Zhang             irn    += count;
1579a6053eceSJunchao Zhang             jcn    += count;
1580a6053eceSJunchao Zhang             val    += count;
1581a6053eceSJunchao Zhang             count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1582a6053eceSJunchao Zhang             remain -= count;
1583a6053eceSJunchao Zhang           }
15843ab56b82SJunchao Zhang         }
15853ab56b82SJunchao Zhang       } else {
1586a6053eceSJunchao Zhang         irn    = mumps->irn;
1587a6053eceSJunchao Zhang         jcn    = mumps->jcn;
1588a6053eceSJunchao Zhang         val    = mumps->val;
1589a6053eceSJunchao Zhang         count  = PetscMin(mumps->nnz,PETSC_MPI_INT_MAX);
1590a6053eceSJunchao Zhang         remain = mumps->nnz - count;
1591a6053eceSJunchao Zhang         while (count>0) {
15929566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(irn,count,MPIU_MUMPSINT,0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
15939566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(jcn,count,MPIU_MUMPSINT,0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
15949566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(val,count,MPIU_SCALAR,  0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
1595a6053eceSJunchao Zhang           irn    += count;
1596a6053eceSJunchao Zhang           jcn    += count;
1597a6053eceSJunchao Zhang           val    += count;
1598a6053eceSJunchao Zhang           count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1599a6053eceSJunchao Zhang           remain -= count;
16003ab56b82SJunchao Zhang         }
16013ab56b82SJunchao Zhang       }
1602a6053eceSJunchao Zhang     } else {
1603a6053eceSJunchao Zhang       nreqs = 0;
1604a6053eceSJunchao Zhang       if (mumps->is_omp_master) {
1605a6053eceSJunchao Zhang         val = mumps->val + mumps->recvcount[0];
1606a6053eceSJunchao Zhang         for (i=1; i<osize; i++) { /* Remote communication only since self data is already in place */
1607a6053eceSJunchao Zhang           count  = PetscMin(mumps->recvcount[i],PETSC_MPI_INT_MAX);
1608a6053eceSJunchao Zhang           remain = mumps->recvcount[i] - count;
1609a6053eceSJunchao Zhang           while (count>0) {
16109566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Irecv(val,count,MPIU_SCALAR,i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
1611a6053eceSJunchao Zhang             val    += count;
1612a6053eceSJunchao Zhang             count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1613a6053eceSJunchao Zhang             remain -= count;
1614a6053eceSJunchao Zhang           }
1615a6053eceSJunchao Zhang         }
1616a6053eceSJunchao Zhang       } else {
1617a6053eceSJunchao Zhang         val    = mumps->val;
1618a6053eceSJunchao Zhang         count  = PetscMin(mumps->nnz,PETSC_MPI_INT_MAX);
1619a6053eceSJunchao Zhang         remain = mumps->nnz - count;
1620a6053eceSJunchao Zhang         while (count>0) {
16219566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(val,count,MPIU_SCALAR,0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]));
1622a6053eceSJunchao Zhang           val    += count;
1623a6053eceSJunchao Zhang           count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1624a6053eceSJunchao Zhang           remain -= count;
1625a6053eceSJunchao Zhang         }
1626a6053eceSJunchao Zhang       }
1627a6053eceSJunchao Zhang     }
16289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nreqs,mumps->reqs,MPI_STATUSES_IGNORE));
1629a6053eceSJunchao Zhang     mumps->tag++; /* It is totally fine for above send/recvs to share one mpi tag */
1630a6053eceSJunchao Zhang   }
16313ab56b82SJunchao Zhang   PetscFunctionReturn(0);
16323ab56b82SJunchao Zhang }
16333ab56b82SJunchao Zhang 
16340481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1635af281ebdSHong Zhang {
1636e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1637ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
1638397b6df1SKris Buschelman 
1639397b6df1SKris Buschelman   PetscFunctionBegin;
1640dbf6bb8dSprj-   if (mumps->id.INFOG(1) < 0 && !(mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0)) {
16412aca8efcSHong Zhang     if (mumps->id.INFOG(1) == -6) {
16429566063dSJacob Faibussowitsch       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)));
16436baea169SHong Zhang     }
16449566063dSJacob 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)));
16452aca8efcSHong Zhang     PetscFunctionReturn(0);
16462aca8efcSHong Zhang   }
16476baea169SHong Zhang 
16489566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, mumps));
16499566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_REUSE_MATRIX,mumps));
1650397b6df1SKris Buschelman 
1651397b6df1SKris Buschelman   /* numerical factorization phase */
1652329ec9b3SHong Zhang   /*-------------------------------*/
1653a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
16544e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1655a5e57a09SHong Zhang     if (!mumps->myid) {
1656940cd9d6SSatish Balay       mumps->id.a = (MumpsScalar*)mumps->val;
1657397b6df1SKris Buschelman     }
1658397b6df1SKris Buschelman   } else {
1659940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1660397b6df1SKris Buschelman   }
16613ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
1662a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
1663c0d63f2fSHong Zhang     if (A->erroriffailure) {
166498921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d",mumps->id.INFOG(1),mumps->id.INFO(2));
1665151787a6SHong Zhang     } else {
1666c0d63f2fSHong Zhang       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
16679566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1668603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1669c0d63f2fSHong Zhang       } else if (mumps->id.INFOG(1) == -13) {
16709566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1671603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1672c0d63f2fSHong Zhang       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10)) {
16739566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1674603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
16752aca8efcSHong Zhang       } else {
16769566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1677603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OTHER;
1678151787a6SHong Zhang       }
16792aca8efcSHong Zhang     }
1680397b6df1SKris Buschelman   }
16812c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mumps->myid && mumps->id.ICNTL(16) > 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"  mumps->id.ICNTL(16):=%d",mumps->id.INFOG(16));
1682397b6df1SKris Buschelman 
1683b3cb21ddSStefano Zampini   F->assembled    = PETSC_TRUE;
1684d47f36abSHong Zhang 
1685b3cb21ddSStefano Zampini   if (F->schur) { /* reset Schur status to unfactored */
16863cb7dd0eSStefano Zampini #if defined(PETSC_HAVE_CUDA)
1687c70f7ee4SJunchao Zhang     F->schur->offloadmask = PETSC_OFFLOAD_CPU;
16883cb7dd0eSStefano Zampini #endif
1689b3cb21ddSStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1690b3cb21ddSStefano Zampini       mumps->id.ICNTL(19) = 2;
16919566063dSJacob Faibussowitsch       PetscCall(MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur));
1692b3cb21ddSStefano Zampini     }
16939566063dSJacob Faibussowitsch     PetscCall(MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED));
1694b3cb21ddSStefano Zampini   }
169567877ebaSShri Abhyankar 
1696066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1697066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1698066565c5SStefano Zampini 
16993ab56b82SJunchao Zhang   if (!mumps->is_omp_master) mumps->id.INFO(23) = 0;
17002d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
170167877ebaSShri Abhyankar     PetscInt    lsol_loc;
170267877ebaSShri Abhyankar     PetscScalar *sol_loc;
17032205254eSKarl Rupp 
17049566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ));
1705c2093ab7SHong Zhang 
1706c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
1707c2093ab7SHong Zhang     if (mumps->x_seq) {
17089566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&mumps->scat_sol));
17099566063dSJacob Faibussowitsch       PetscCall(PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc));
17109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mumps->x_seq));
1711c2093ab7SHong Zhang     }
1712a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
17139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc));
1714a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
1715940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
17169566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq));
171767877ebaSShri Abhyankar   }
17189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(mumps->id.RINFO(2)));
1719397b6df1SKris Buschelman   PetscFunctionReturn(0);
1720397b6df1SKris Buschelman }
1721397b6df1SKris Buschelman 
17229a2535b5SHong Zhang /* Sets MUMPS options from the options database */
17239a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1724dcd589f8SShri Abhyankar {
1725e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1726dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1727a6053eceSJunchao Zhang   PetscMUMPSInt  icntl=0;
1728a6053eceSJunchao Zhang   PetscInt       info[80],i,ninfo=80;
1729a6053eceSJunchao Zhang   PetscBool      flg=PETSC_FALSE;
1730dcd589f8SShri Abhyankar 
1731dcd589f8SShri Abhyankar   PetscFunctionBegin;
17329566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");PetscCall(ierr);
17339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg));
17349a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
17359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg));
17369a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
17379566063dSJacob 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));
17389a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
1739dcd589f8SShri Abhyankar 
17409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg));
17419a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
17429a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
17439a2535b5SHong Zhang 
17449566063dSJacob 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));
17459a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
17469a2535b5SHong Zhang 
17479566063dSJacob 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));
1748dcd589f8SShri Abhyankar   if (flg) {
17492c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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");
1750b53c1a7fSBarry Smith     mumps->id.ICNTL(7) = icntl;
1751dcd589f8SShri Abhyankar   }
1752e0b74bf9SHong Zhang 
17539566063dSJacob 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));
17549566063dSJacob 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() */
17559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL));
17569566063dSJacob 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));
17579566063dSJacob 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));
17589566063dSJacob 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));
17599566063dSJacob 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));
17609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL));
176159ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
17629566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F->schur));
17639566063dSJacob Faibussowitsch     PetscCall(MatMumpsResetSchur_Private(mumps));
176459ac8732SStefano Zampini   }
176525aac85cSJunchao Zhang 
176643f3b051SJunchao Zhang   /* Two MPICH Fortran MPI_IN_PLACE binding bugs prevented the use of 'mpich + mumps'. One happened with "mpi4py + mpich + mumps",
176743f3b051SJunchao Zhang      and was reported by Firedrake. See https://bitbucket.org/mpi4py/mpi4py/issues/162/mpi4py-initialization-breaks-fortran
176825aac85cSJunchao Zhang      and a petsc-maint mailing list thread with subject 'MUMPS segfaults in parallel because of ...'
176943f3b051SJunchao Zhang      This bug was fixed by https://github.com/pmodels/mpich/pull/4149. But the fix brought a new bug,
177043f3b051SJunchao Zhang      see https://github.com/pmodels/mpich/issues/5589. This bug was fixed by https://github.com/pmodels/mpich/pull/5590.
177143f3b051SJunchao Zhang      In short, we could not use distributed RHS with MPICH until v4.0b1.
177225aac85cSJunchao Zhang    */
177343f3b051SJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_LT(5,3,0) || (defined(PETSC_HAVE_MPICH_NUMVERSION) && (PETSC_HAVE_MPICH_NUMVERSION < 40000101))
177425aac85cSJunchao Zhang   mumps->ICNTL20 = 0;  /* Centralized dense RHS*/
177543f3b051SJunchao Zhang #else
177643f3b051SJunchao Zhang   mumps->ICNTL20 = 10; /* Distributed dense RHS*/
177725aac85cSJunchao Zhang #endif
17789566063dSJacob 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));
17792c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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);
178025aac85cSJunchao Zhang #if PETSC_PKG_MUMPS_VERSION_LT(5,3,0)
17812c71b3e2SJacob Faibussowitsch   PetscCheckFalse(flg && mumps->ICNTL20 == 10,PETSC_COMM_SELF,PETSC_ERR_SUP,"ICNTL(20)=10 is not supported before MUMPS-5.3.0");
178225aac85cSJunchao Zhang #endif
17839566063dSJacob 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 */
17849a2535b5SHong Zhang 
17859566063dSJacob 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));
17869566063dSJacob 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));
17879566063dSJacob 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));
17889a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
17899a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1790d7ebd59bSHong Zhang   }
1791d7ebd59bSHong Zhang 
17929566063dSJacob 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));
17939566063dSJacob 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));
17949566063dSJacob 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));
17959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL));
17969566063dSJacob 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));
17979566063dSJacob 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 */
17989566063dSJacob 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));
17999566063dSJacob Faibussowitsch   /* PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL));  -- not supported by PETSc API */
18009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL));
18019566063dSJacob 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));
18029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsMUMPSInt("-mat_mumps_icntl_36","ICNTL(36): choice of BLR factorization variant","None",mumps->id.ICNTL(36),&mumps->id.ICNTL(36),NULL));
18039566063dSJacob 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));
1804dcd589f8SShri Abhyankar 
18059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL));
18069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL));
18079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL));
18089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL));
18099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL));
18109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL));
1811e5bb22a1SHong Zhang 
18129566063dSJacob 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));
1813b34f08ffSHong Zhang 
18149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL));
1815b34f08ffSHong Zhang   if (ninfo) {
1816*08401ef6SPierre Jolivet     PetscCheck(ninfo <= 80,PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %" PetscInt_FMT " must <= 80",ninfo);
18179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ninfo,&mumps->info));
1818b34f08ffSHong Zhang     mumps->ninfo = ninfo;
1819b34f08ffSHong Zhang     for (i=0; i<ninfo; i++) {
18202c71b3e2SJacob Faibussowitsch       PetscCheckFalse(info[i] < 0 || info[i]>80,PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %" PetscInt_FMT " must between 1 and 80",ninfo);
18212a808120SBarry Smith       else  mumps->info[i] = info[i];
1822b34f08ffSHong Zhang     }
1823b34f08ffSHong Zhang   }
1824b34f08ffSHong Zhang 
18259566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
1826dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1827dcd589f8SShri Abhyankar }
1828dcd589f8SShri Abhyankar 
1829f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1830dcd589f8SShri Abhyankar {
18317c405c4aSJunchao Zhang   PetscInt       nthreads=0;
1832dcd589f8SShri Abhyankar 
1833dcd589f8SShri Abhyankar   PetscFunctionBegin;
18343ab56b82SJunchao Zhang   mumps->petsc_comm = PetscObjectComm((PetscObject)A);
18359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(mumps->petsc_comm,&mumps->petsc_size));
18369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(mumps->petsc_comm,&mumps->myid));/* "if (!myid)" still works even if mumps_comm is different */
18373ab56b82SJunchao Zhang 
18389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,((PetscObject)A)->prefix,"-mat_mumps_use_omp_threads",&mumps->use_petsc_omp_support));
18397c405c4aSJunchao Zhang   if (mumps->use_petsc_omp_support) nthreads = -1; /* -1 will let PetscOmpCtrlCreate() guess a proper value when user did not supply one */
18409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,((PetscObject)A)->prefix,"-mat_mumps_use_omp_threads",&nthreads,NULL));
18413ab56b82SJunchao Zhang   if (mumps->use_petsc_omp_support) {
18423ab56b82SJunchao Zhang #if defined(PETSC_HAVE_OPENMP_SUPPORT)
18439566063dSJacob Faibussowitsch     PetscCall(PetscOmpCtrlCreate(mumps->petsc_comm,nthreads,&mumps->omp_ctrl));
18449566063dSJacob Faibussowitsch     PetscCall(PetscOmpCtrlGetOmpComms(mumps->omp_ctrl,&mumps->omp_comm,&mumps->mumps_comm,&mumps->is_omp_master));
18453ab56b82SJunchao Zhang #else
184698921bdaSJacob Faibussowitsch     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",((PetscObject)A)->prefix?((PetscObject)A)->prefix:"");
18473ab56b82SJunchao Zhang #endif
18483ab56b82SJunchao Zhang   } else {
18493ab56b82SJunchao Zhang     mumps->omp_comm      = PETSC_COMM_SELF;
18503ab56b82SJunchao Zhang     mumps->mumps_comm    = mumps->petsc_comm;
18513ab56b82SJunchao Zhang     mumps->is_omp_master = PETSC_TRUE;
18523ab56b82SJunchao Zhang   }
18539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(mumps->omp_comm,&mumps->omp_comm_size));
1854a6053eceSJunchao Zhang   mumps->reqs = NULL;
1855a6053eceSJunchao Zhang   mumps->tag  = 0;
18562205254eSKarl Rupp 
185725aac85cSJunchao Zhang   /* It looks like MUMPS does not dup the input comm. Dup a new comm for MUMPS to avoid any tag mismatches. */
185825aac85cSJunchao Zhang   if (mumps->mumps_comm != MPI_COMM_NULL) {
18599566063dSJacob Faibussowitsch     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)A),&mumps->mumps_comm));
186025aac85cSJunchao Zhang   }
186125aac85cSJunchao Zhang 
18622d4298aeSJunchao Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->mumps_comm);
1863f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
1864f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1865f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
18663ab56b82SJunchao Zhang 
18673ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
1868*08401ef6SPierre Jolivet   PetscCheck(mumps->id.INFOG(1) >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in PetscInitializeMUMPS: INFOG(1)=%d",mumps->id.INFOG(1));
18693ab56b82SJunchao Zhang 
18703ab56b82SJunchao Zhang   /* copy MUMPS default control values from master to slaves. Although slaves do not call MUMPS, they may access these values in code.
18713ab56b82SJunchao Zhang      For example, ICNTL(9) is initialized to 1 by MUMPS and slaves check ICNTL(9) in MatSolve_MUMPS.
18723ab56b82SJunchao Zhang    */
18739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(mumps->id.icntl,40,MPI_INT,  0,mumps->omp_comm));
18749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(mumps->id.cntl, 15,MPIU_REAL,0,mumps->omp_comm));
1875f697e70eSHong Zhang 
18760298fd71SBarry Smith   mumps->scat_rhs = NULL;
18770298fd71SBarry Smith   mumps->scat_sol = NULL;
18789a2535b5SHong Zhang 
187970544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
18809a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
18819a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
18822d4298aeSJunchao Zhang   if (mumps->petsc_size == 1) {
18839a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1884ff8cf4dfSJunchao Zhang     mumps->id.ICNTL(7)  = 7;   /* automatic choice of ordering done by the package */
18859a2535b5SHong Zhang   } else {
18869a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
188770544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
18889a2535b5SHong Zhang   }
18896444a565SStefano Zampini 
18906444a565SStefano Zampini   /* schur */
18916444a565SStefano Zampini   mumps->id.size_schur    = 0;
18926444a565SStefano Zampini   mumps->id.listvar_schur = NULL;
18936444a565SStefano Zampini   mumps->id.schur         = NULL;
1894b5fa320bSStefano Zampini   mumps->sizeredrhs       = 0;
189559ac8732SStefano Zampini   mumps->schur_sol        = NULL;
189659ac8732SStefano Zampini   mumps->schur_sizesol    = 0;
1897dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1898dcd589f8SShri Abhyankar }
1899dcd589f8SShri Abhyankar 
19009a625307SHong Zhang PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
19015cd7cf9dSHong Zhang {
19025cd7cf9dSHong Zhang   PetscFunctionBegin;
19035cd7cf9dSHong Zhang   if (mumps->id.INFOG(1) < 0) {
19045cd7cf9dSHong Zhang     if (A->erroriffailure) {
190598921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d",mumps->id.INFOG(1));
19065cd7cf9dSHong Zhang     } else {
19075cd7cf9dSHong Zhang       if (mumps->id.INFOG(1) == -6) {
19089566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1909603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
19105cd7cf9dSHong Zhang       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
19119566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1912603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1913dbf6bb8dSprj-       } else if (mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0) {
19149566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"Empty matrix\n"));
19155cd7cf9dSHong Zhang       } else {
19169566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)));
1917603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OTHER;
19185cd7cf9dSHong Zhang       }
19195cd7cf9dSHong Zhang     }
19205cd7cf9dSHong Zhang   }
19215cd7cf9dSHong Zhang   PetscFunctionReturn(0);
19225cd7cf9dSHong Zhang }
19235cd7cf9dSHong Zhang 
19240481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1925b24902e0SBarry Smith {
1926e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
192767877ebaSShri Abhyankar   Vec            b;
192867877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1929397b6df1SKris Buschelman 
1930397b6df1SKris Buschelman   PetscFunctionBegin;
1931d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
1932d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
1933d47f36abSHong Zhang     PetscFunctionReturn(0);
1934d47f36abSHong Zhang   }
1935dcd589f8SShri Abhyankar 
19369a2535b5SHong Zhang   /* Set MUMPS options from the options database */
19379566063dSJacob Faibussowitsch   PetscCall(PetscSetMUMPSFromOptions(F,A));
1938dcd589f8SShri Abhyankar 
19399566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
19409566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX,mumps));
1941dcd589f8SShri Abhyankar 
194267877ebaSShri Abhyankar   /* analysis phase */
194367877ebaSShri Abhyankar   /*----------------*/
1944a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1945a5e57a09SHong Zhang   mumps->id.n   = M;
1946a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
194767877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1948a5e57a09SHong Zhang     if (!mumps->myid) {
1949a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
1950a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
1951a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
1952a6053eceSJunchao Zhang       if (mumps->id.ICNTL(6)>1) mumps->id.a = (MumpsScalar*)mumps->val;
19534ac6704cSBarry Smith       if (r) {
19544ac6704cSBarry Smith         mumps->id.ICNTL(7) = 1;
1955a5e57a09SHong Zhang         if (!mumps->myid) {
1956e0b74bf9SHong Zhang           const PetscInt *idx;
1957a6053eceSJunchao Zhang           PetscInt       i;
19582205254eSKarl Rupp 
19599566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(M,&mumps->id.perm_in));
19609566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(r,&idx));
19619566063dSJacob Faibussowitsch           for (i=0; i<M; i++) PetscCall(PetscMUMPSIntCast(idx[i]+1,&(mumps->id.perm_in[i]))); /* perm_in[]: start from 1, not 0! */
19629566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(r,&idx));
1963e0b74bf9SHong Zhang         }
1964e0b74bf9SHong Zhang       }
196567877ebaSShri Abhyankar     }
196667877ebaSShri Abhyankar     break;
196767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1968a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
1969a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
1970a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
1971a6053eceSJunchao Zhang     if (mumps->id.ICNTL(6)>1) mumps->id.a_loc = (MumpsScalar*)mumps->val;
197225aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
19739566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A,NULL,&b));
19749566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq));
19759566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
197625aac85cSJunchao Zhang     }
197767877ebaSShri Abhyankar     break;
197867877ebaSShri Abhyankar   }
19793ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
19809566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps));
198167877ebaSShri Abhyankar 
1982719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1983dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
198451d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
19854e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1986eb3ef3b2SHong Zhang   F->ops->mattransposesolve = MatMatTransposeSolve_MUMPS;
1987d47f36abSHong Zhang 
1988d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
1989b24902e0SBarry Smith   PetscFunctionReturn(0);
1990b24902e0SBarry Smith }
1991b24902e0SBarry Smith 
1992450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1993450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1994450b117fSShri Abhyankar {
1995e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
199667877ebaSShri Abhyankar   Vec            b;
199767877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1998450b117fSShri Abhyankar 
1999450b117fSShri Abhyankar   PetscFunctionBegin;
2000d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2001d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
2002d47f36abSHong Zhang     PetscFunctionReturn(0);
2003d47f36abSHong Zhang   }
2004dcd589f8SShri Abhyankar 
20059a2535b5SHong Zhang   /* Set MUMPS options from the options database */
20069566063dSJacob Faibussowitsch   PetscCall(PetscSetMUMPSFromOptions(F,A));
2007dcd589f8SShri Abhyankar 
20089566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
20099566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX,mumps));
201067877ebaSShri Abhyankar 
201167877ebaSShri Abhyankar   /* analysis phase */
201267877ebaSShri Abhyankar   /*----------------*/
2013a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
2014a5e57a09SHong Zhang   mumps->id.n   = M;
2015a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
201667877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
2017a5e57a09SHong Zhang     if (!mumps->myid) {
2018a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2019a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2020a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2021a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
2022940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
202367877ebaSShri Abhyankar       }
202467877ebaSShri Abhyankar     }
202567877ebaSShri Abhyankar     break;
202667877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
2027a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2028a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2029a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2030a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
2031940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
203267877ebaSShri Abhyankar     }
203325aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
20349566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A,NULL,&b));
20359566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq));
20369566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
203725aac85cSJunchao Zhang     }
203867877ebaSShri Abhyankar     break;
203967877ebaSShri Abhyankar   }
20403ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
20419566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps));
204267877ebaSShri Abhyankar 
2043450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
2044dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
204551d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
2046d47f36abSHong Zhang 
2047d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
2048450b117fSShri Abhyankar   PetscFunctionReturn(0);
2049450b117fSShri Abhyankar }
2050b24902e0SBarry Smith 
2051141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
205267877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
2053b24902e0SBarry Smith {
2054e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
205567877ebaSShri Abhyankar   Vec            b;
205667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
2057397b6df1SKris Buschelman 
2058397b6df1SKris Buschelman   PetscFunctionBegin;
2059d47f36abSHong Zhang   if (mumps->matstruc == SAME_NONZERO_PATTERN) {
2060d47f36abSHong Zhang     /* F is assembled by a previous call of MatLUFactorSymbolic_AIJMUMPS() */
2061d47f36abSHong Zhang     PetscFunctionReturn(0);
2062d47f36abSHong Zhang   }
2063dcd589f8SShri Abhyankar 
20649a2535b5SHong Zhang   /* Set MUMPS options from the options database */
20659566063dSJacob Faibussowitsch   PetscCall(PetscSetMUMPSFromOptions(F,A));
2066dcd589f8SShri Abhyankar 
20679566063dSJacob Faibussowitsch   PetscCall((*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps));
20689566063dSJacob Faibussowitsch   PetscCall(MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX,mumps));
2069dcd589f8SShri Abhyankar 
207067877ebaSShri Abhyankar   /* analysis phase */
207167877ebaSShri Abhyankar   /*----------------*/
2072a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
2073a5e57a09SHong Zhang   mumps->id.n   = M;
2074a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
207567877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
2076a5e57a09SHong Zhang     if (!mumps->myid) {
2077a6053eceSJunchao Zhang       mumps->id.nnz = mumps->nnz;
2078a6053eceSJunchao Zhang       mumps->id.irn = mumps->irn;
2079a6053eceSJunchao Zhang       mumps->id.jcn = mumps->jcn;
2080a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
2081940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
208267877ebaSShri Abhyankar       }
208367877ebaSShri Abhyankar     }
208467877ebaSShri Abhyankar     break;
208567877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
2086a6053eceSJunchao Zhang     mumps->id.nnz_loc = mumps->nnz;
2087a6053eceSJunchao Zhang     mumps->id.irn_loc = mumps->irn;
2088a6053eceSJunchao Zhang     mumps->id.jcn_loc = mumps->jcn;
2089a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
2090940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
209167877ebaSShri Abhyankar     }
209225aac85cSJunchao Zhang     if (mumps->ICNTL20 == 0) { /* Centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
20939566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A,NULL,&b));
20949566063dSJacob Faibussowitsch       PetscCall(VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq));
20959566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
209625aac85cSJunchao Zhang     }
209767877ebaSShri Abhyankar     break;
209867877ebaSShri Abhyankar   }
20993ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
21009566063dSJacob Faibussowitsch   PetscCall(MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps));
21015cd7cf9dSHong Zhang 
21022792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
2103dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
210451d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
21054e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
210623a5080aSHong Zhang   F->ops->mattransposesolve     = MatMatTransposeSolve_MUMPS;
21074e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
21080298fd71SBarry Smith   F->ops->getinertia = NULL;
21094e34a73bSHong Zhang #else
21104e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
2111db4efbfdSBarry Smith #endif
2112d47f36abSHong Zhang 
2113d47f36abSHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
2114b24902e0SBarry Smith   PetscFunctionReturn(0);
2115b24902e0SBarry Smith }
2116b24902e0SBarry Smith 
211764e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
211874ed9c26SBarry Smith {
211964e6c443SBarry Smith   PetscBool         iascii;
212064e6c443SBarry Smith   PetscViewerFormat format;
2121e69c285eSBarry Smith   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
2122f6c57405SHong Zhang 
2123f6c57405SHong Zhang   PetscFunctionBegin;
212464e6c443SBarry Smith   /* check if matrix is mumps type */
212564e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
212664e6c443SBarry Smith 
21279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
212864e6c443SBarry Smith   if (iascii) {
21299566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
213064e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
21319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n"));
21329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d\n",mumps->id.sym));
21339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d\n",mumps->id.par));
21349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d\n",mumps->id.ICNTL(1)));
21359566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d\n",mumps->id.ICNTL(2)));
21369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d\n",mumps->id.ICNTL(3)));
21379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d\n",mumps->id.ICNTL(4)));
21389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d\n",mumps->id.ICNTL(5)));
21399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d\n",mumps->id.ICNTL(6)));
21409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d\n",mumps->id.ICNTL(7)));
21419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d\n",mumps->id.ICNTL(8)));
21429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d\n",mumps->id.ICNTL(10)));
21439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d\n",mumps->id.ICNTL(11)));
2144a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
21459566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4)));
21469566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5)));
21479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6)));
21489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8)));
21499566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9)));
21509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11)));
2151f6c57405SHong Zhang       }
21529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d\n",mumps->id.ICNTL(12)));
21539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (sequential factorization of the root node):  %d\n",mumps->id.ICNTL(13)));
21549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d\n",mumps->id.ICNTL(14)));
2155f6c57405SHong Zhang       /* ICNTL(15-17) not used */
21569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d\n",mumps->id.ICNTL(18)));
21579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                      %d\n",mumps->id.ICNTL(19)));
21589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (RHS sparse pattern):                         %d\n",mumps->id.ICNTL(20)));
21599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d\n",mumps->id.ICNTL(21)));
21609566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d\n",mumps->id.ICNTL(22)));
21619566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d\n",mumps->id.ICNTL(23)));
2162c0165424SHong Zhang 
21639566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d\n",mumps->id.ICNTL(24)));
21649566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d\n",mumps->id.ICNTL(25)));
21659566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for RHS or solution):          %d\n",mumps->id.ICNTL(26)));
21669566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (blocking size for multiple RHS):             %d\n",mumps->id.ICNTL(27)));
21679566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d\n",mumps->id.ICNTL(28)));
21689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d\n",mumps->id.ICNTL(29)));
216942179a6aSHong Zhang 
21709566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d\n",mumps->id.ICNTL(30)));
21719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d\n",mumps->id.ICNTL(31)));
21729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d\n",mumps->id.ICNTL(33)));
21739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(35) (activate BLR based factorization):           %d\n",mumps->id.ICNTL(35)));
21749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(36) (choice of BLR factorization variant):        %d\n",mumps->id.ICNTL(36)));
21759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  ICNTL(38) (estimated compression rate of LU factors):   %d\n",mumps->id.ICNTL(38)));
2176f6c57405SHong Zhang 
21779566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1)));
21789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2)));
21799566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3)));
21809566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4)));
21819566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5)));
21829566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7)));
2183f6c57405SHong Zhang 
2184a5b23f4aSJose E. Roman       /* information local to each processor */
21859566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n"));
21869566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
21879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1)));
21889566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
21899566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n"));
21909566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2)));
21919566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
21929566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n"));
21939566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3)));
21949566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2195f6c57405SHong Zhang 
21969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n"));
21979566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d\n",mumps->myid,mumps->id.INFO(15)));
21989566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2199f6c57405SHong Zhang 
22009566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n"));
22019566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d\n",mumps->myid,mumps->id.INFO(16)));
22029566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2203f6c57405SHong Zhang 
22049566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n"));
22059566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d\n",mumps->myid,mumps->id.INFO(23)));
22069566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2207b34f08ffSHong Zhang 
2208a0e18203SThibaut Appel       if (mumps->ninfo && mumps->ninfo <= 80) {
2209b34f08ffSHong Zhang         PetscInt i;
2210b34f08ffSHong Zhang         for (i=0; i<mumps->ninfo; i++) {
22119566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  INFO(%" PetscInt_FMT "): \n",mumps->info[i]));
22129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d\n",mumps->myid,mumps->id.INFO(mumps->info[i])));
22139566063dSJacob Faibussowitsch           PetscCall(PetscViewerFlush(viewer));
2214b34f08ffSHong Zhang         }
2215b34f08ffSHong Zhang       }
22169566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2217f6c57405SHong Zhang 
2218a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
22199566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1)));
22209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2)));
22219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3)));
22229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34)));
2223f6c57405SHong Zhang 
22249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d\n",mumps->id.INFOG(3)));
22259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d\n",mumps->id.INFOG(4)));
22269566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d\n",mumps->id.INFOG(5)));
22279566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d\n",mumps->id.INFOG(6)));
22289566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively used after analysis): %d\n",mumps->id.INFOG(7)));
22299566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d\n",mumps->id.INFOG(8)));
22309566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d\n",mumps->id.INFOG(9)));
22319566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d\n",mumps->id.INFOG(10)));
22329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d\n",mumps->id.INFOG(11)));
22339566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d\n",mumps->id.INFOG(12)));
22349566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d\n",mumps->id.INFOG(13)));
22359566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d\n",mumps->id.INFOG(14)));
22369566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d\n",mumps->id.INFOG(15)));
22379566063dSJacob 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)));
22389566063dSJacob 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)));
22399566063dSJacob 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)));
22409566063dSJacob 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)));
22419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d\n",mumps->id.INFOG(20)));
22429566063dSJacob 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)));
22439566063dSJacob 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)));
22449566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d\n",mumps->id.INFOG(23)));
22459566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d\n",mumps->id.INFOG(24)));
22469566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d\n",mumps->id.INFOG(25)));
22479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28)));
22489566063dSJacob 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)));
22499566063dSJacob 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)));
22509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32)));
22519566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33)));
22529566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34)));
22539566063dSJacob 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)));
22549566063dSJacob 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)));
22559566063dSJacob 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)));
22569566063dSJacob 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)));
22579566063dSJacob 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)));
2258f6c57405SHong Zhang       }
2259f6c57405SHong Zhang     }
2260cb828f0fSHong Zhang   }
2261f6c57405SHong Zhang   PetscFunctionReturn(0);
2262f6c57405SHong Zhang }
2263f6c57405SHong Zhang 
226435bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
226535bd34faSBarry Smith {
2266e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;
226735bd34faSBarry Smith 
226835bd34faSBarry Smith   PetscFunctionBegin;
226935bd34faSBarry Smith   info->block_size        = 1.0;
2270cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
2271cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
227235bd34faSBarry Smith   info->nz_unneeded       = 0.0;
227335bd34faSBarry Smith   info->assemblies        = 0.0;
227435bd34faSBarry Smith   info->mallocs           = 0.0;
227535bd34faSBarry Smith   info->memory            = 0.0;
227635bd34faSBarry Smith   info->fill_ratio_given  = 0;
227735bd34faSBarry Smith   info->fill_ratio_needed = 0;
227835bd34faSBarry Smith   info->factor_mallocs    = 0;
227935bd34faSBarry Smith   PetscFunctionReturn(0);
228035bd34faSBarry Smith }
228135bd34faSBarry Smith 
22825ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
22838e7ba810SStefano Zampini PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
22846444a565SStefano Zampini {
2285e69c285eSBarry Smith   Mat_MUMPS         *mumps =(Mat_MUMPS*)F->data;
2286a3d589ffSStefano Zampini   const PetscScalar *arr;
22878e7ba810SStefano Zampini   const PetscInt    *idxs;
22888e7ba810SStefano Zampini   PetscInt          size,i;
22896444a565SStefano Zampini 
22906444a565SStefano Zampini   PetscFunctionBegin;
22919566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&size));
22922d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
22933ab56b82SJunchao Zhang     PetscBool ls,gs; /* gs is false if any rank other than root has non-empty IS */
2294241dbb5eSStefano Zampini 
22953ab56b82SJunchao Zhang     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE; /* always true on root; false on others if their size != 0 */
22969566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->petsc_comm));
229728b400f6SJacob Faibussowitsch     PetscCheck(gs,PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc");
2298241dbb5eSStefano Zampini   }
2299b3cb21ddSStefano Zampini 
2300b3cb21ddSStefano Zampini   /* Schur complement matrix */
23019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F->schur));
23029566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,size,size,NULL,&F->schur));
23039566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(F->schur,&arr));
2304a3d589ffSStefano Zampini   mumps->id.schur      = (MumpsScalar*)arr;
2305a3d589ffSStefano Zampini   mumps->id.size_schur = size;
2306a3d589ffSStefano Zampini   mumps->id.schur_lld  = size;
23079566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(F->schur,&arr));
2308b3cb21ddSStefano Zampini   if (mumps->sym == 1) {
23099566063dSJacob Faibussowitsch     PetscCall(MatSetOption(F->schur,MAT_SPD,PETSC_TRUE));
2310b3cb21ddSStefano Zampini   }
2311b3cb21ddSStefano Zampini 
2312b3cb21ddSStefano Zampini   /* MUMPS expects Fortran style indices */
23139566063dSJacob Faibussowitsch   PetscCall(PetscFree(mumps->id.listvar_schur));
23149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&mumps->id.listvar_schur));
23159566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is,&idxs));
23169566063dSJacob Faibussowitsch   for (i=0; i<size; i++) PetscCall(PetscMUMPSIntCast(idxs[i]+1,&(mumps->id.listvar_schur[i])));
23179566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is,&idxs));
23182d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
2319241dbb5eSStefano Zampini     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
2320241dbb5eSStefano Zampini   } else {
23216444a565SStefano Zampini     if (F->factortype == MAT_FACTOR_LU) {
232259ac8732SStefano Zampini       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
23236444a565SStefano Zampini     } else {
232459ac8732SStefano Zampini       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
23256444a565SStefano Zampini     }
2326241dbb5eSStefano Zampini   }
232759ac8732SStefano Zampini   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
2328b5fa320bSStefano Zampini   mumps->id.ICNTL(26) = -1;
23296444a565SStefano Zampini   PetscFunctionReturn(0);
23306444a565SStefano Zampini }
233159ac8732SStefano Zampini 
23326444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
23335a05ddb0SStefano Zampini PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
23346444a565SStefano Zampini {
23356444a565SStefano Zampini   Mat            St;
2336e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
23376444a565SStefano Zampini   PetscScalar    *array;
23386444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX)
23398ac429a0SStefano Zampini   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
23406444a565SStefano Zampini #endif
23416444a565SStefano Zampini 
23426444a565SStefano Zampini   PetscFunctionBegin;
2343*08401ef6SPierre Jolivet   PetscCheck(mumps->id.ICNTL(19),PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
23449566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&St));
23459566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur));
23469566063dSJacob Faibussowitsch   PetscCall(MatSetType(St,MATDENSE));
23479566063dSJacob Faibussowitsch   PetscCall(MatSetUp(St));
23489566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(St,&array));
234959ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full matrix */
23506444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
23516444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
23526444a565SStefano Zampini       for (i=0;i<N;i++) {
23536444a565SStefano Zampini         for (j=0;j<N;j++) {
23546444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
23556444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
23566444a565SStefano Zampini #else
23576444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
23586444a565SStefano Zampini #endif
23596444a565SStefano Zampini           array[j*N+i] = val;
23606444a565SStefano Zampini         }
23616444a565SStefano Zampini       }
23626444a565SStefano Zampini     } else { /* stored by columns */
23639566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur));
23646444a565SStefano Zampini     }
23656444a565SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
23666444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
23676444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
23686444a565SStefano Zampini       for (i=0;i<N;i++) {
23696444a565SStefano Zampini         for (j=i;j<N;j++) {
23706444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
23716444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
23726444a565SStefano Zampini #else
23736444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
23746444a565SStefano Zampini #endif
23756444a565SStefano Zampini           array[i*N+j] = val;
23766444a565SStefano Zampini           array[j*N+i] = val;
23776444a565SStefano Zampini         }
23786444a565SStefano Zampini       }
23796444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
23809566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur));
23816444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
23826444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
23836444a565SStefano Zampini       for (i=0;i<N;i++) {
23846444a565SStefano Zampini         for (j=0;j<i+1;j++) {
23856444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
23866444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
23876444a565SStefano Zampini #else
23886444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
23896444a565SStefano Zampini #endif
23906444a565SStefano Zampini           array[i*N+j] = val;
23916444a565SStefano Zampini           array[j*N+i] = val;
23926444a565SStefano Zampini         }
23936444a565SStefano Zampini       }
23946444a565SStefano Zampini     }
23956444a565SStefano Zampini   }
23969566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(St,&array));
23976444a565SStefano Zampini   *S   = St;
23986444a565SStefano Zampini   PetscFunctionReturn(0);
23996444a565SStefano Zampini }
24006444a565SStefano Zampini 
240159ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
24025ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
24035ccb76cbSHong Zhang {
2404e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
24055ccb76cbSHong Zhang 
24065ccb76cbSHong Zhang   PetscFunctionBegin;
24079566063dSJacob Faibussowitsch   PetscCall(PetscMUMPSIntCast(ival,&mumps->id.ICNTL(icntl)));
24085ccb76cbSHong Zhang   PetscFunctionReturn(0);
24095ccb76cbSHong Zhang }
24105ccb76cbSHong Zhang 
2411bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2412bc6112feSHong Zhang {
2413e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2414bc6112feSHong Zhang 
2415bc6112feSHong Zhang   PetscFunctionBegin;
2416bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
2417bc6112feSHong Zhang   PetscFunctionReturn(0);
2418bc6112feSHong Zhang }
2419bc6112feSHong Zhang 
24205ccb76cbSHong Zhang /*@
24215ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
24225ccb76cbSHong Zhang 
24235ccb76cbSHong Zhang    Logically Collective on Mat
24245ccb76cbSHong Zhang 
24255ccb76cbSHong Zhang    Input Parameters:
24265ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
24275ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
24285ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
24295ccb76cbSHong Zhang 
24305ccb76cbSHong Zhang   Options Database:
2431147403d9SBarry Smith .   -mat_mumps_icntl_<icntl> <ival> - change the option numbered icntl to ival
24325ccb76cbSHong Zhang 
24335ccb76cbSHong Zhang    Level: beginner
24345ccb76cbSHong Zhang 
243596a0c994SBarry Smith    References:
2436606c0280SSatish Balay .  * - MUMPS Users' Guide
24375ccb76cbSHong Zhang 
243822235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
24395ccb76cbSHong Zhang @*/
24405ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
24415ccb76cbSHong Zhang {
24425ccb76cbSHong Zhang   PetscFunctionBegin;
24432989dfd4SHong Zhang   PetscValidType(F,1);
244428b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
24455ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
24465ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
2447cac4c232SBarry Smith   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
24485ccb76cbSHong Zhang   PetscFunctionReturn(0);
24495ccb76cbSHong Zhang }
24505ccb76cbSHong Zhang 
2451a21f80fcSHong Zhang /*@
2452a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
2453a21f80fcSHong Zhang 
2454a21f80fcSHong Zhang    Logically Collective on Mat
2455a21f80fcSHong Zhang 
2456a21f80fcSHong Zhang    Input Parameters:
2457a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2458a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
2459a21f80fcSHong Zhang 
2460a21f80fcSHong Zhang   Output Parameter:
2461a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
2462a21f80fcSHong Zhang 
2463a21f80fcSHong Zhang    Level: beginner
2464a21f80fcSHong Zhang 
246596a0c994SBarry Smith    References:
2466606c0280SSatish Balay .  * - MUMPS Users' Guide
2467a21f80fcSHong Zhang 
246822235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2469a21f80fcSHong Zhang @*/
2470bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2471bc6112feSHong Zhang {
2472bc6112feSHong Zhang   PetscFunctionBegin;
24732989dfd4SHong Zhang   PetscValidType(F,1);
247428b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2475bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2476bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
2477cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2478bc6112feSHong Zhang   PetscFunctionReturn(0);
2479bc6112feSHong Zhang }
2480bc6112feSHong Zhang 
24818928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
24828928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
24838928b65cSHong Zhang {
2484e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
24858928b65cSHong Zhang 
24868928b65cSHong Zhang   PetscFunctionBegin;
24878928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
24888928b65cSHong Zhang   PetscFunctionReturn(0);
24898928b65cSHong Zhang }
24908928b65cSHong Zhang 
2491bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2492bc6112feSHong Zhang {
2493e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2494bc6112feSHong Zhang 
2495bc6112feSHong Zhang   PetscFunctionBegin;
2496bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
2497bc6112feSHong Zhang   PetscFunctionReturn(0);
2498bc6112feSHong Zhang }
2499bc6112feSHong Zhang 
25008928b65cSHong Zhang /*@
25018928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
25028928b65cSHong Zhang 
25038928b65cSHong Zhang    Logically Collective on Mat
25048928b65cSHong Zhang 
25058928b65cSHong Zhang    Input Parameters:
25068928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
25078928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
25088928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
25098928b65cSHong Zhang 
25108928b65cSHong Zhang   Options Database:
2511147403d9SBarry Smith .   -mat_mumps_cntl_<icntl> <val>  - change the option numbered icntl to ival
25128928b65cSHong Zhang 
25138928b65cSHong Zhang    Level: beginner
25148928b65cSHong Zhang 
251596a0c994SBarry Smith    References:
2516606c0280SSatish Balay .  * - MUMPS Users' Guide
25178928b65cSHong Zhang 
251822235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
25198928b65cSHong Zhang @*/
25208928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
25218928b65cSHong Zhang {
25228928b65cSHong Zhang   PetscFunctionBegin;
25232989dfd4SHong Zhang   PetscValidType(F,1);
252428b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
25258928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2526bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
2527cac4c232SBarry Smith   PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
25288928b65cSHong Zhang   PetscFunctionReturn(0);
25298928b65cSHong Zhang }
25308928b65cSHong Zhang 
2531a21f80fcSHong Zhang /*@
2532a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
2533a21f80fcSHong Zhang 
2534a21f80fcSHong Zhang    Logically Collective on Mat
2535a21f80fcSHong Zhang 
2536a21f80fcSHong Zhang    Input Parameters:
2537a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2538a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
2539a21f80fcSHong Zhang 
2540a21f80fcSHong Zhang   Output Parameter:
2541a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
2542a21f80fcSHong Zhang 
2543a21f80fcSHong Zhang    Level: beginner
2544a21f80fcSHong Zhang 
254596a0c994SBarry Smith    References:
2546606c0280SSatish Balay .  * - MUMPS Users' Guide
2547a21f80fcSHong Zhang 
254822235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2549a21f80fcSHong Zhang @*/
2550bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2551bc6112feSHong Zhang {
2552bc6112feSHong Zhang   PetscFunctionBegin;
25532989dfd4SHong Zhang   PetscValidType(F,1);
255428b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2555bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2556bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2557cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2558bc6112feSHong Zhang   PetscFunctionReturn(0);
2559bc6112feSHong Zhang }
2560bc6112feSHong Zhang 
2561ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2562bc6112feSHong Zhang {
2563e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2564bc6112feSHong Zhang 
2565bc6112feSHong Zhang   PetscFunctionBegin;
2566bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
2567bc6112feSHong Zhang   PetscFunctionReturn(0);
2568bc6112feSHong Zhang }
2569bc6112feSHong Zhang 
2570ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2571bc6112feSHong Zhang {
2572e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2573bc6112feSHong Zhang 
2574bc6112feSHong Zhang   PetscFunctionBegin;
2575bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
2576bc6112feSHong Zhang   PetscFunctionReturn(0);
2577bc6112feSHong Zhang }
2578bc6112feSHong Zhang 
2579ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2580bc6112feSHong Zhang {
2581e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2582bc6112feSHong Zhang 
2583bc6112feSHong Zhang   PetscFunctionBegin;
2584bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
2585bc6112feSHong Zhang   PetscFunctionReturn(0);
2586bc6112feSHong Zhang }
2587bc6112feSHong Zhang 
2588ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2589bc6112feSHong Zhang {
2590e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2591bc6112feSHong Zhang 
2592bc6112feSHong Zhang   PetscFunctionBegin;
2593bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
2594bc6112feSHong Zhang   PetscFunctionReturn(0);
2595bc6112feSHong Zhang }
2596bc6112feSHong Zhang 
259789a9c03aSHong Zhang PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS)
2598bb599dfdSHong Zhang {
25990e6b8875SHong Zhang   Mat            Bt = NULL,Btseq = NULL;
26000e6b8875SHong Zhang   PetscBool      flg;
2601bb599dfdSHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
2602bb599dfdSHong Zhang   PetscScalar    *aa;
2603f410b75aSHong Zhang   PetscInt       spnr,*ia,*ja,M,nrhs;
2604bb599dfdSHong Zhang 
2605bb599dfdSHong Zhang   PetscFunctionBegin;
2606064a246eSJacob Faibussowitsch   PetscValidPointer(spRHS,2);
26079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flg));
26080e6b8875SHong Zhang   if (flg) {
26099566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(spRHS,&Bt));
26100e6b8875SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix");
2611bb599dfdSHong Zhang 
26129566063dSJacob Faibussowitsch   PetscCall(MatMumpsSetIcntl(F,30,1));
2613bb599dfdSHong Zhang 
26142d4298aeSJunchao Zhang   if (mumps->petsc_size > 1) {
26150e6b8875SHong Zhang     Mat_MPIAIJ *b = (Mat_MPIAIJ*)Bt->data;
26160e6b8875SHong Zhang     Btseq = b->A;
26170e6b8875SHong Zhang   } else {
26180e6b8875SHong Zhang     Btseq = Bt;
26190e6b8875SHong Zhang   }
26200e6b8875SHong Zhang 
26219566063dSJacob Faibussowitsch   PetscCall(MatGetSize(spRHS,&M,&nrhs));
2622f410b75aSHong Zhang   mumps->id.nrhs = nrhs;
2623f410b75aSHong Zhang   mumps->id.lrhs = M;
2624f410b75aSHong Zhang   mumps->id.rhs  = NULL;
2625f410b75aSHong Zhang 
2626e3f2db6aSHong Zhang   if (!mumps->myid) {
26279566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArray(Btseq,&aa));
26289566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(Btseq,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg));
262928b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
26309566063dSJacob Faibussowitsch     PetscCall(PetscMUMPSIntCSRCast(mumps,spnr,ia,ja,&mumps->id.irhs_ptr,&mumps->id.irhs_sparse,&mumps->id.nz_rhs));
2631bb599dfdSHong Zhang     mumps->id.rhs_sparse  = (MumpsScalar*)aa;
2632e3f2db6aSHong Zhang   } else {
2633e3f2db6aSHong Zhang     mumps->id.irhs_ptr    = NULL;
2634e3f2db6aSHong Zhang     mumps->id.irhs_sparse = NULL;
2635e3f2db6aSHong Zhang     mumps->id.nz_rhs      = 0;
2636e3f2db6aSHong Zhang     mumps->id.rhs_sparse  = NULL;
2637e3f2db6aSHong Zhang   }
2638bb599dfdSHong Zhang   mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
2639e3f2db6aSHong Zhang   mumps->id.ICNTL(21)   = 0; /* solution is in assembled centralized format */
2640bb599dfdSHong Zhang 
2641bb599dfdSHong Zhang   /* solve phase */
2642bb599dfdSHong Zhang   /*-------------*/
2643bb599dfdSHong Zhang   mumps->id.job = JOB_SOLVE;
26443ab56b82SJunchao Zhang   PetscMUMPS_c(mumps);
2645e3f2db6aSHong Zhang   if (mumps->id.INFOG(1) < 0)
264698921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d INFO(2)=%d",mumps->id.INFOG(1),mumps->id.INFO(2));
264714267174SHong Zhang 
2648e3f2db6aSHong Zhang   if (!mumps->myid) {
26499566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJRestoreArray(Btseq,&aa));
26509566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(Btseq,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg));
265128b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
2652e3f2db6aSHong Zhang   }
2653bb599dfdSHong Zhang   PetscFunctionReturn(0);
2654bb599dfdSHong Zhang }
2655bb599dfdSHong Zhang 
2656bb599dfdSHong Zhang /*@
265789a9c03aSHong Zhang   MatMumpsGetInverse - Get user-specified set of entries in inverse of A
2658bb599dfdSHong Zhang 
2659bb599dfdSHong Zhang    Logically Collective on Mat
2660bb599dfdSHong Zhang 
2661bb599dfdSHong Zhang    Input Parameters:
2662bb599dfdSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2663e3f2db6aSHong Zhang -  spRHS - sequential sparse matrix in MATTRANSPOSEMAT format holding specified indices in processor[0]
2664bb599dfdSHong Zhang 
2665bb599dfdSHong Zhang   Output Parameter:
2666e3f2db6aSHong Zhang . spRHS - requested entries of inverse of A
2667bb599dfdSHong Zhang 
2668bb599dfdSHong Zhang    Level: beginner
2669bb599dfdSHong Zhang 
2670bb599dfdSHong Zhang    References:
2671606c0280SSatish Balay .  * - MUMPS Users' Guide
2672bb599dfdSHong Zhang 
2673bb599dfdSHong Zhang .seealso: MatGetFactor(), MatCreateTranspose()
2674bb599dfdSHong Zhang @*/
267589a9c03aSHong Zhang PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS)
2676bb599dfdSHong Zhang {
2677bb599dfdSHong Zhang   PetscFunctionBegin;
2678bb599dfdSHong Zhang   PetscValidType(F,1);
267928b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2680cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));
2681bb599dfdSHong Zhang   PetscFunctionReturn(0);
2682bb599dfdSHong Zhang }
2683bb599dfdSHong Zhang 
26840e6b8875SHong Zhang PetscErrorCode MatMumpsGetInverseTranspose_MUMPS(Mat F,Mat spRHST)
26850e6b8875SHong Zhang {
26860e6b8875SHong Zhang   Mat            spRHS;
26870e6b8875SHong Zhang 
26880e6b8875SHong Zhang   PetscFunctionBegin;
26899566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(spRHST,&spRHS));
26909566063dSJacob Faibussowitsch   PetscCall(MatMumpsGetInverse_MUMPS(F,spRHS));
26919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&spRHS));
26920e6b8875SHong Zhang   PetscFunctionReturn(0);
26930e6b8875SHong Zhang }
26940e6b8875SHong Zhang 
26950e6b8875SHong Zhang /*@
2696eef1237cSHong Zhang   MatMumpsGetInverseTranspose - Get user-specified set of entries in inverse of matrix A^T
26970e6b8875SHong Zhang 
26980e6b8875SHong Zhang    Logically Collective on Mat
26990e6b8875SHong Zhang 
27000e6b8875SHong Zhang    Input Parameters:
27010e6b8875SHong Zhang +  F - the factored matrix of A obtained by calling MatGetFactor() from PETSc-MUMPS interface
27020e6b8875SHong Zhang -  spRHST - sequential sparse matrix in MATAIJ format holding specified indices of A^T in processor[0]
27030e6b8875SHong Zhang 
27040e6b8875SHong Zhang   Output Parameter:
27050e6b8875SHong Zhang . spRHST - requested entries of inverse of A^T
27060e6b8875SHong Zhang 
27070e6b8875SHong Zhang    Level: beginner
27080e6b8875SHong Zhang 
27090e6b8875SHong Zhang    References:
2710606c0280SSatish Balay .  * - MUMPS Users' Guide
27110e6b8875SHong Zhang 
27120e6b8875SHong Zhang .seealso: MatGetFactor(), MatCreateTranspose(), MatMumpsGetInverse()
27130e6b8875SHong Zhang @*/
27140e6b8875SHong Zhang PetscErrorCode MatMumpsGetInverseTranspose(Mat F,Mat spRHST)
27150e6b8875SHong Zhang {
27160e6b8875SHong Zhang   PetscBool      flg;
27170e6b8875SHong Zhang 
27180e6b8875SHong Zhang   PetscFunctionBegin;
27190e6b8875SHong Zhang   PetscValidType(F,1);
272028b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
27219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)spRHST,&flg,MATSEQAIJ,MATMPIAIJ,NULL));
272228b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)spRHST),PETSC_ERR_ARG_WRONG,"Matrix spRHST must be MATAIJ matrix");
27230e6b8875SHong Zhang 
2724cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetInverseTranspose_C",(Mat,Mat),(F,spRHST));
27250e6b8875SHong Zhang   PetscFunctionReturn(0);
27260e6b8875SHong Zhang }
27270e6b8875SHong Zhang 
2728a21f80fcSHong Zhang /*@
2729a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
2730a21f80fcSHong Zhang 
2731a21f80fcSHong Zhang    Logically Collective on Mat
2732a21f80fcSHong Zhang 
2733a21f80fcSHong Zhang    Input Parameters:
2734a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2735a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
2736a21f80fcSHong Zhang 
2737a21f80fcSHong Zhang   Output Parameter:
2738a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
2739a21f80fcSHong Zhang 
2740a21f80fcSHong Zhang    Level: beginner
2741a21f80fcSHong Zhang 
274296a0c994SBarry Smith    References:
2743606c0280SSatish Balay .  * - MUMPS Users' Guide
2744a21f80fcSHong Zhang 
274522235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2746a21f80fcSHong Zhang @*/
2747ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2748bc6112feSHong Zhang {
2749bc6112feSHong Zhang   PetscFunctionBegin;
27502989dfd4SHong Zhang   PetscValidType(F,1);
275128b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2752ca810319SHong Zhang   PetscValidIntPointer(ival,3);
2753cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2754bc6112feSHong Zhang   PetscFunctionReturn(0);
2755bc6112feSHong Zhang }
2756bc6112feSHong Zhang 
2757a21f80fcSHong Zhang /*@
2758a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2759a21f80fcSHong Zhang 
2760a21f80fcSHong Zhang    Logically Collective on Mat
2761a21f80fcSHong Zhang 
2762a21f80fcSHong Zhang    Input Parameters:
2763a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2764a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
2765a21f80fcSHong Zhang 
2766a21f80fcSHong Zhang   Output Parameter:
2767a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
2768a21f80fcSHong Zhang 
2769a21f80fcSHong Zhang    Level: beginner
2770a21f80fcSHong Zhang 
277196a0c994SBarry Smith    References:
2772606c0280SSatish Balay .  * - MUMPS Users' Guide
2773a21f80fcSHong Zhang 
277422235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2775a21f80fcSHong Zhang @*/
2776ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2777bc6112feSHong Zhang {
2778bc6112feSHong Zhang   PetscFunctionBegin;
27792989dfd4SHong Zhang   PetscValidType(F,1);
278028b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2781ca810319SHong Zhang   PetscValidIntPointer(ival,3);
2782cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2783bc6112feSHong Zhang   PetscFunctionReturn(0);
2784bc6112feSHong Zhang }
2785bc6112feSHong Zhang 
2786a21f80fcSHong Zhang /*@
2787a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2788a21f80fcSHong Zhang 
2789a21f80fcSHong Zhang    Logically Collective on Mat
2790a21f80fcSHong Zhang 
2791a21f80fcSHong Zhang    Input Parameters:
2792a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2793a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
2794a21f80fcSHong Zhang 
2795a21f80fcSHong Zhang   Output Parameter:
2796a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
2797a21f80fcSHong Zhang 
2798a21f80fcSHong Zhang    Level: beginner
2799a21f80fcSHong Zhang 
280096a0c994SBarry Smith    References:
2801606c0280SSatish Balay .  * - MUMPS Users' Guide
2802a21f80fcSHong Zhang 
280322235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfog()
2804a21f80fcSHong Zhang @*/
2805ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2806bc6112feSHong Zhang {
2807bc6112feSHong Zhang   PetscFunctionBegin;
28082989dfd4SHong Zhang   PetscValidType(F,1);
280928b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2810bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2811cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2812bc6112feSHong Zhang   PetscFunctionReturn(0);
2813bc6112feSHong Zhang }
2814bc6112feSHong Zhang 
2815a21f80fcSHong Zhang /*@
2816a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2817a21f80fcSHong Zhang 
2818a21f80fcSHong Zhang    Logically Collective on Mat
2819a21f80fcSHong Zhang 
2820a21f80fcSHong Zhang    Input Parameters:
2821a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2822a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
2823a21f80fcSHong Zhang 
2824a21f80fcSHong Zhang   Output Parameter:
2825a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
2826a21f80fcSHong Zhang 
2827a21f80fcSHong Zhang    Level: beginner
2828a21f80fcSHong Zhang 
282996a0c994SBarry Smith    References:
2830606c0280SSatish Balay .  * - MUMPS Users' Guide
2831a21f80fcSHong Zhang 
283222235d61SPierre Jolivet .seealso: MatGetFactor(), MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo()
2833a21f80fcSHong Zhang @*/
2834ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2835bc6112feSHong Zhang {
2836bc6112feSHong Zhang   PetscFunctionBegin;
28372989dfd4SHong Zhang   PetscValidType(F,1);
283828b400f6SJacob Faibussowitsch   PetscCheck(F->factortype,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2839bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2840cac4c232SBarry Smith   PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2841bc6112feSHong Zhang   PetscFunctionReturn(0);
2842bc6112feSHong Zhang }
2843bc6112feSHong Zhang 
284424b6179bSKris Buschelman /*MC
28452692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
284624b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
284724b6179bSKris Buschelman 
284841c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
284924b6179bSKris Buschelman 
2850c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS
2851c2b89b5dSBarry Smith 
2852217d3b1eSJunchao Zhang   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. See details below.
2853217d3b1eSJunchao Zhang 
28543ca39a21SBarry Smith   Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver
2855c2b89b5dSBarry Smith 
285624b6179bSKris Buschelman   Options Database Keys:
28574422a9fcSPatrick Sanan +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
28584422a9fcSPatrick Sanan .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
28594422a9fcSPatrick Sanan .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
28604422a9fcSPatrick Sanan .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
28614422a9fcSPatrick Sanan .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2862b53c1a7fSBarry 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
2863b53c1a7fSBarry Smith                         Use -pc_factor_mat_ordering_type <type> to have PETSc perform the ordering (sequential only)
28644422a9fcSPatrick Sanan .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
28654422a9fcSPatrick Sanan .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
28664422a9fcSPatrick Sanan .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
28674422a9fcSPatrick Sanan .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
28684422a9fcSPatrick Sanan .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
28694422a9fcSPatrick Sanan .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
28704422a9fcSPatrick Sanan .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
287125aac85cSJunchao Zhang .  -mat_mumps_icntl_20  - ICNTL(20): give MUMPS centralized (0) or distributed (10) dense RHS
28724422a9fcSPatrick Sanan .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
28734422a9fcSPatrick Sanan .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
28744422a9fcSPatrick Sanan .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
28754422a9fcSPatrick Sanan .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
28764422a9fcSPatrick Sanan .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
28774422a9fcSPatrick Sanan .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
28784422a9fcSPatrick Sanan .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
28794422a9fcSPatrick Sanan .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
28804422a9fcSPatrick Sanan .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
28814422a9fcSPatrick Sanan .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2882a0e18203SThibaut Appel .  -mat_mumps_icntl_35 - ICNTL(35): level of activation of BLR (Block Low-Rank) feature
2883a0e18203SThibaut Appel .  -mat_mumps_icntl_36 - ICNTL(36): controls the choice of BLR factorization variant
2884a0e18203SThibaut Appel .  -mat_mumps_icntl_38 - ICNTL(38): sets the estimated compression rate of LU factors with BLR
28854422a9fcSPatrick Sanan .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
28864422a9fcSPatrick Sanan .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
28874422a9fcSPatrick Sanan .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
28884422a9fcSPatrick Sanan .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2889217d3b1eSJunchao Zhang .  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2890a0e18203SThibaut Appel .  -mat_mumps_cntl_7 - CNTL(7): precision of the dropping parameter used during BLR factorization
2891217d3b1eSJunchao 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.
2892217d3b1eSJunchao Zhang                                    Default might be the number of cores per CPU package (socket) as reported by hwloc and suggested by the MUMPS manual.
289324b6179bSKris Buschelman 
289424b6179bSKris Buschelman   Level: beginner
289524b6179bSKris Buschelman 
289695452b02SPatrick Sanan     Notes:
289738548759SBarry Smith     MUMPS Cholesky does not handle (complex) Hermitian matrices http://mumps.enseeiht.fr/doc/userguide_5.2.1.pdf so using it will error if the matrix is Hermitian.
289838548759SBarry Smith 
2899c0decd05SBarry 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 the failure by calling
29009fc87aa7SBarry Smith $          KSPGetPC(ksp,&pc);
29019fc87aa7SBarry Smith $          PCFactorGetMatrix(pc,&mat);
29029fc87aa7SBarry Smith $          MatMumpsGetInfo(mat,....);
29039fc87aa7SBarry Smith $          MatMumpsGetInfog(mat,....); etc.
29049fc87aa7SBarry Smith            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.
29059fc87aa7SBarry Smith 
2906a5399872SJunchao Zhang   Using MUMPS with 64-bit integers
2907a5399872SJunchao Zhang     MUMPS provides 64-bit integer support in two build modes:
2908a5399872SJunchao Zhang       full 64-bit: here MUMPS is built with C preprocessing flag -DINTSIZE64 and Fortran compiler option -i8, -fdefault-integer-8 or equivalent, and
2909a5399872SJunchao 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).
29108fcaa860SBarry Smith 
2911a5399872SJunchao 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,
2912a5399872SJunchao 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
2913a5399872SJunchao 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
2914a5399872SJunchao Zhang       integer) build of all dependent libraries MPI, ScaLAPACK, LAPACK and BLAS.
2915a5399872SJunchao Zhang 
2916a5399872SJunchao 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.
2917a5399872SJunchao Zhang 
2918a5399872SJunchao Zhang   Two modes to run MUMPS/PETSc with OpenMP
29198fcaa860SBarry Smith $     Set OMP_NUM_THREADS and run with fewer MPI ranks than cores. For example, if you want to have 16 OpenMP
29208fcaa860SBarry Smith $     threads per rank, then you may use "export OMP_NUM_THREADS=16 && mpirun -n 4 ./test".
29218fcaa860SBarry Smith 
29228fcaa860SBarry Smith $     -mat_mumps_use_omp_threads [m] and run your code with as many MPI ranks as the number of cores. For example,
29238fcaa860SBarry 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"
29248fcaa860SBarry Smith 
29258fcaa860SBarry Smith    To run MUMPS in MPI+OpenMP hybrid mode (i.e., enable multithreading in MUMPS), but still run the non-MUMPS part
2926217d3b1eSJunchao Zhang    (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
2927217d3b1eSJunchao Zhang    (or --with-hwloc), and have an MPI that supports MPI-3.0's process shared memory (which is usually available). Since MUMPS calls BLAS
29288fcaa860SBarry Smith    libraries, to really get performance, you should have multithreaded BLAS libraries such as Intel MKL, AMD ACML, Cray libSci or OpenBLAS
29298fcaa860SBarry Smith    (PETSc will automatically try to utilized a threaded BLAS if --with-openmp is provided).
2930217d3b1eSJunchao Zhang 
29318fcaa860SBarry 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
2932217d3b1eSJunchao Zhang    processes on each compute node. Listing the processes in rank ascending order, we split processes on a node into consecutive groups of
2933217d3b1eSJunchao 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
2934217d3b1eSJunchao 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
2935217d3b1eSJunchao 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.
2936217d3b1eSJunchao 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,
2937217d3b1eSJunchao 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
2938217d3b1eSJunchao Zhang    MPI ranks to cores, then with -mat_mumps_use_omp_threads 16, a master rank (and threads it spawns) will use half cores in socket 0, and half
2939217d3b1eSJunchao 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
2940217d3b1eSJunchao Zhang    problem will not happen. Therefore, when you use -mat_mumps_use_omp_threads, you need to keep an eye on your MPI rank mapping and CPU binding.
29418fcaa860SBarry 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
2942217d3b1eSJunchao Zhang    examine the mapping result.
2943217d3b1eSJunchao Zhang 
2944217d3b1eSJunchao Zhang    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,
2945217d3b1eSJunchao Zhang    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
2946217d3b1eSJunchao Zhang    calls omp_set_num_threads(m) internally before calling MUMPS.
2947217d3b1eSJunchao Zhang 
2948217d3b1eSJunchao Zhang    References:
2949606c0280SSatish Balay +  * - Heroux, Michael A., R. Brightwell, and Michael M. Wolf. "Bi-modal MPI and MPI+ threads computing on scalable multicore systems." IJHPCA (Submitted) (2011).
2950606c0280SSatish Balay -  * - Gutierrez, Samuel K., et al. "Accommodating Thread-Level Heterogeneity in Coupled Parallel Applications." Parallel and Distributed Processing Symposium (IPDPS), 2017 IEEE International. IEEE, 2017.
2951217d3b1eSJunchao Zhang 
2952e24fcbf7SPierre Jolivet .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetIcntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCFactorGetMatrix()
295341c8de11SBarry Smith 
295424b6179bSKris Buschelman M*/
295524b6179bSKris Buschelman 
2956ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
295735bd34faSBarry Smith {
295835bd34faSBarry Smith   PetscFunctionBegin;
29592692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
296035bd34faSBarry Smith   PetscFunctionReturn(0);
296135bd34faSBarry Smith }
296235bd34faSBarry Smith 
2963bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
2964cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
29652877fffaSHong Zhang {
29662877fffaSHong Zhang   Mat            B;
29672877fffaSHong Zhang   Mat_MUMPS      *mumps;
2968ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
29692c7c0729SBarry Smith   PetscMPIInt    size;
29702877fffaSHong Zhang 
29712877fffaSHong Zhang   PetscFunctionBegin;
2972eb1ec7c1SStefano Zampini  #if defined(PETSC_USE_COMPLEX)
29732c71b3e2SJacob Faibussowitsch   PetscCheckFalse(A->hermitian && !A->symmetric && ftype == MAT_FACTOR_CHOLESKY,PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian CHOLESKY Factor is not supported");
2974eb1ec7c1SStefano Zampini  #endif
29752877fffaSHong Zhang   /* Create the factorization matrix */
29769566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ));
29779566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
29789566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
29799566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps",&((PetscObject)B)->type_name));
29809566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
29812877fffaSHong Zhang 
29829566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&mumps));
29832205254eSKarl Rupp 
29842877fffaSHong Zhang   B->ops->view    = MatView_MUMPS;
298535bd34faSBarry Smith   B->ops->getinfo = MatGetInfo_MUMPS;
29862205254eSKarl Rupp 
29879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps));
29889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS));
29899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS));
29909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS));
29919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS));
29929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS));
29939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS));
29949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS));
29959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS));
29969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS));
29979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS));
29989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS));
29999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS));
30006444a565SStefano Zampini 
3001450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
3002450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
3003d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
3004bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
3005bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
30069566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]));
3007746480a1SHong Zhang     mumps->sym = 0;
3008dcd589f8SShri Abhyankar   } else {
300967877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
3010450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
3011bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
3012bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
30139566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
301459ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
301559ac8732SStefano Zampini     mumps->sym = 2;
301659ac8732SStefano Zampini #else
30176fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
30186fdc2a6dSBarry Smith     else                      mumps->sym = 2;
301959ac8732SStefano Zampini #endif
3020450b117fSShri Abhyankar   }
30212877fffaSHong Zhang 
302200c67f3bSHong Zhang   /* set solvertype */
30239566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
30249566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype));
30259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
30262c7c0729SBarry Smith   if (size == 1) {
30274ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3028f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
30292c7c0729SBarry Smith   }
30302877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
3031e69c285eSBarry Smith   B->data         = (void*)mumps;
30322205254eSKarl Rupp 
30339566063dSJacob Faibussowitsch   PetscCall(PetscInitializeMUMPS(A,mumps));
3034746480a1SHong Zhang 
30352877fffaSHong Zhang   *F = B;
3036d47f36abSHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
30372877fffaSHong Zhang   PetscFunctionReturn(0);
30382877fffaSHong Zhang }
30392877fffaSHong Zhang 
3040bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
3041cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
30422877fffaSHong Zhang {
30432877fffaSHong Zhang   Mat            B;
30442877fffaSHong Zhang   Mat_MUMPS      *mumps;
3045ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
30462c7c0729SBarry Smith   PetscMPIInt    size;
30472877fffaSHong Zhang 
30482877fffaSHong Zhang   PetscFunctionBegin;
3049eb1ec7c1SStefano Zampini  #if defined(PETSC_USE_COMPLEX)
3050*08401ef6SPierre Jolivet   PetscCheck(!A->hermitian || A->symmetric,PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian CHOLESKY Factor is not supported");
3051eb1ec7c1SStefano Zampini  #endif
30529566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
30539566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
30549566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps",&((PetscObject)B)->type_name));
30559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3056e69c285eSBarry Smith 
30579566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&mumps));
30589566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ));
3059bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
306016ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
3061dcd589f8SShri Abhyankar   } else {
3062bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
3063bccb9932SShri Abhyankar   }
3064bccb9932SShri Abhyankar 
306567877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
3066bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
3067722b6324SPierre Jolivet   B->ops->getinfo                = MatGetInfo_MUMPS;
30682205254eSKarl Rupp 
30699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps));
30709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS));
30719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS));
30729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS));
30739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS));
30749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS));
30759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS));
30769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS));
30779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS));
30789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS));
30799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS));
30809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS));
30819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS));
30822205254eSKarl Rupp 
3083f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
308459ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
308559ac8732SStefano Zampini   mumps->sym = 2;
308659ac8732SStefano Zampini #else
30876fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
30886fdc2a6dSBarry Smith   else                      mumps->sym = 2;
308959ac8732SStefano Zampini #endif
3090a214ac2aSShri Abhyankar 
309100c67f3bSHong Zhang   /* set solvertype */
30929566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
30939566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype));
30949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
30952c7c0729SBarry Smith   if (size == 1) {
30964ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3097f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
30982c7c0729SBarry Smith   }
30999566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
3100f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
3101e69c285eSBarry Smith   B->data         = (void*)mumps;
31022205254eSKarl Rupp 
31039566063dSJacob Faibussowitsch   PetscCall(PetscInitializeMUMPS(A,mumps));
3104746480a1SHong Zhang 
31052877fffaSHong Zhang   *F = B;
3106d47f36abSHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
31072877fffaSHong Zhang   PetscFunctionReturn(0);
31082877fffaSHong Zhang }
310997969023SHong Zhang 
3110cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
311167877ebaSShri Abhyankar {
311267877ebaSShri Abhyankar   Mat            B;
311367877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
3114ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
31152c7c0729SBarry Smith   PetscMPIInt    size;
311667877ebaSShri Abhyankar 
311767877ebaSShri Abhyankar   PetscFunctionBegin;
311867877ebaSShri Abhyankar   /* Create the factorization matrix */
31199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ));
31209566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
31219566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
31229566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps",&((PetscObject)B)->type_name));
31239566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3124450b117fSShri Abhyankar 
31259566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&mumps));
3126450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
3127450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
3128450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
3129bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
3130bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
3131746480a1SHong Zhang     mumps->sym = 0;
31329566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]));
3133546078acSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead");
3134bccb9932SShri Abhyankar 
3135450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
3136722b6324SPierre Jolivet   B->ops->getinfo     = MatGetInfo_MUMPS;
31372205254eSKarl Rupp 
31389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps));
31399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS));
31409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS));
31419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS));
31429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS));
31439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS));
31449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS));
31459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS));
31469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS));
31479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS));
31489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS));
31499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS));
31509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS));
3151450b117fSShri Abhyankar 
315200c67f3bSHong Zhang   /* set solvertype */
31539566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
31549566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype));
31559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
31562c7c0729SBarry Smith   if (size == 1) {
31574ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization */
3158f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
31592c7c0729SBarry Smith   }
31607ee00b23SStefano Zampini   B->ops->destroy = MatDestroy_MUMPS;
31617ee00b23SStefano Zampini   B->data         = (void*)mumps;
31627ee00b23SStefano Zampini 
31639566063dSJacob Faibussowitsch   PetscCall(PetscInitializeMUMPS(A,mumps));
31647ee00b23SStefano Zampini 
31657ee00b23SStefano Zampini   *F = B;
3166d47f36abSHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
31677ee00b23SStefano Zampini   PetscFunctionReturn(0);
31687ee00b23SStefano Zampini }
31697ee00b23SStefano Zampini 
31707ee00b23SStefano Zampini /* MatGetFactor for Seq and MPI SELL matrices */
31717ee00b23SStefano Zampini static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F)
31727ee00b23SStefano Zampini {
31737ee00b23SStefano Zampini   Mat            B;
31747ee00b23SStefano Zampini   Mat_MUMPS      *mumps;
31757ee00b23SStefano Zampini   PetscBool      isSeqSELL;
31762c7c0729SBarry Smith   PetscMPIInt    size;
31777ee00b23SStefano Zampini 
31787ee00b23SStefano Zampini   PetscFunctionBegin;
31797ee00b23SStefano Zampini   /* Create the factorization matrix */
31809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL));
31819566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
31829566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
31839566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("mumps",&((PetscObject)B)->type_name));
31849566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
31857ee00b23SStefano Zampini 
31869566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&mumps));
31877ee00b23SStefano Zampini 
31887ee00b23SStefano Zampini   B->ops->view        = MatView_MUMPS;
31897ee00b23SStefano Zampini   B->ops->getinfo     = MatGetInfo_MUMPS;
31907ee00b23SStefano Zampini 
31919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps));
31929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS));
31939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS));
31949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS));
31959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS));
31969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS));
31979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS));
31989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS));
31999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS));
32009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS));
32019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS));
32027ee00b23SStefano Zampini 
32037ee00b23SStefano Zampini   if (ftype == MAT_FACTOR_LU) {
32047ee00b23SStefano Zampini     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
32057ee00b23SStefano Zampini     B->factortype            = MAT_FACTOR_LU;
32067ee00b23SStefano Zampini     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
32077ee00b23SStefano Zampini     else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
32087ee00b23SStefano Zampini     mumps->sym = 0;
32099566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]));
32107ee00b23SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
32117ee00b23SStefano Zampini 
32127ee00b23SStefano Zampini   /* set solvertype */
32139566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
32149566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype));
32159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
32162c7c0729SBarry Smith   if (size == 1) {
32174ac6704cSBarry Smith     /* MUMPS option -mat_mumps_icntl_7 1 is automatically set if PETSc ordering is passed into symbolic factorization  */
3218f73b0415SBarry Smith     B->canuseordering = PETSC_TRUE;
32192c7c0729SBarry Smith   }
3220450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
3221e69c285eSBarry Smith   B->data         = (void*)mumps;
32222205254eSKarl Rupp 
32239566063dSJacob Faibussowitsch   PetscCall(PetscInitializeMUMPS(A,mumps));
3224746480a1SHong Zhang 
3225450b117fSShri Abhyankar   *F = B;
3226d47f36abSHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
3227450b117fSShri Abhyankar   PetscFunctionReturn(0);
3228450b117fSShri Abhyankar }
322942c9c57cSBarry Smith 
32303ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
323142c9c57cSBarry Smith {
323242c9c57cSBarry Smith   PetscFunctionBegin;
32339566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps));
32349566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps));
32359566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps));
32369566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps));
32379566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps));
32389566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps));
32399566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps));
32409566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps));
32419566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps));
32429566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps));
32439566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps));
324442c9c57cSBarry Smith   PetscFunctionReturn(0);
324542c9c57cSBarry Smith }
3246