xref: /petsc/src/ksp/pc/impls/mpi/pcmpi.c (revision 4568237604e340b5e35202a3c58b4896aff61a77)
1f1f2ae84SBarry Smith /*
2f1f2ae84SBarry Smith     This file creates an MPI parallel KSP from a sequential PC that lives on MPI rank 0.
3f1f2ae84SBarry Smith     It is intended to allow using PETSc MPI parallel linear solvers from non-MPI codes.
4f1f2ae84SBarry Smith 
5dd8e379bSPierre Jolivet     That program may use OpenMP to compute the right-hand side and matrix for the linear system
6f1f2ae84SBarry Smith 
7f1f2ae84SBarry Smith     The code uses MPI_COMM_WORLD below but maybe it should be PETSC_COMM_WORLD
8f1f2ae84SBarry Smith 
9f1f2ae84SBarry Smith     The resulting KSP and PC can only be controlled via the options database, though some common commands
10f1f2ae84SBarry Smith     could be passed through the server.
11f1f2ae84SBarry Smith 
12f1f2ae84SBarry Smith */
13f1f2ae84SBarry Smith #include <petsc/private/pcimpl.h>
14f1f2ae84SBarry Smith #include <petsc/private/kspimpl.h>
15789afff4SPierre Jolivet #include <petscts.h>
16789afff4SPierre Jolivet #include <petsctao.h>
17f1f2ae84SBarry Smith 
18f1f2ae84SBarry Smith #define PC_MPI_MAX_RANKS  256
19f1f2ae84SBarry Smith #define PC_MPI_COMM_WORLD MPI_COMM_WORLD
20f1f2ae84SBarry Smith 
21f1f2ae84SBarry Smith typedef struct {
22f1f2ae84SBarry Smith   KSP         ksps[PC_MPI_MAX_RANKS];                               /* The addresses of the MPI parallel KSP on each rank, NULL when not on a rank. */
23f1f2ae84SBarry Smith   PetscMPIInt sendcount[PC_MPI_MAX_RANKS], displ[PC_MPI_MAX_RANKS]; /* For scatter/gather of rhs/solution */
24f1f2ae84SBarry Smith   PetscMPIInt NZ[PC_MPI_MAX_RANKS], NZdispl[PC_MPI_MAX_RANKS];      /* For scatter of nonzero values in matrix (and nonzero column indices initially */
25f1f2ae84SBarry Smith   PetscInt    mincntperrank;                                        /* minimum number of desired nonzeros per active rank in MPI parallel KSP solve */
26f1f2ae84SBarry Smith   PetscBool   alwaysuseserver;                                      /* for debugging use the server infrastructure even if only one MPI rank is used for the solve */
27f1f2ae84SBarry Smith } PC_MPI;
28f1f2ae84SBarry Smith 
299371c9d4SSatish Balay typedef enum {
309371c9d4SSatish Balay   PCMPI_EXIT, /* exit the PC server loop, means the controlling sequential program is done */
31f1f2ae84SBarry Smith   PCMPI_CREATE,
32f1f2ae84SBarry Smith   PCMPI_SET_MAT,           /* set original matrix (or one with different nonzero pattern) */
33f1f2ae84SBarry Smith   PCMPI_UPDATE_MAT_VALUES, /* update current matrix with new nonzero values */
34f1f2ae84SBarry Smith   PCMPI_SOLVE,
35f1f2ae84SBarry Smith   PCMPI_VIEW,
36f1f2ae84SBarry Smith   PCMPI_DESTROY /* destroy a KSP that is no longer needed */
37f1f2ae84SBarry Smith } PCMPICommand;
38f1f2ae84SBarry Smith 
39f1f2ae84SBarry Smith static MPI_Comm      PCMPIComms[PC_MPI_MAX_RANKS];
40f1f2ae84SBarry Smith static PetscBool     PCMPICommSet = PETSC_FALSE;
41f1f2ae84SBarry Smith static PetscInt      PCMPISolveCounts[PC_MPI_MAX_RANKS], PCMPIKSPCounts[PC_MPI_MAX_RANKS], PCMPIMatCounts[PC_MPI_MAX_RANKS], PCMPISolveCountsSeq = 0, PCMPIKSPCountsSeq = 0;
425316cbedSBarry Smith static PetscInt      PCMPIIterations[PC_MPI_MAX_RANKS], PCMPISizes[PC_MPI_MAX_RANKS], PCMPIIterationsSeq = 0, PCMPISizesSeq = 0;
43*45682376SBarry Smith static PetscLogEvent EventServerDist;
44f1f2ae84SBarry Smith 
45d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPICommsCreate(void)
46d71ae5a4SJacob Faibussowitsch {
47f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
48f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
49f1f2ae84SBarry Smith 
50f1f2ae84SBarry Smith   PetscFunctionBegin;
51f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
52f1f2ae84SBarry Smith   PetscCheck(size <= PC_MPI_MAX_RANKS, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for using more than PC_MPI_MAX_RANKS MPI ranks in an MPI linear solver server solve");
53f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
54f1f2ae84SBarry Smith   /* comm for size 1 is useful only for debugging */
55f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
56f1f2ae84SBarry Smith     PetscMPIInt color = rank < i + 1 ? 0 : MPI_UNDEFINED;
57f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_split(comm, color, 0, &PCMPIComms[i]));
58f1f2ae84SBarry Smith     PCMPISolveCounts[i] = 0;
59f1f2ae84SBarry Smith     PCMPIKSPCounts[i]   = 0;
605316cbedSBarry Smith     PCMPIIterations[i]  = 0;
615316cbedSBarry Smith     PCMPISizes[i]       = 0;
62f1f2ae84SBarry Smith   }
63f1f2ae84SBarry Smith   PCMPICommSet = PETSC_TRUE;
643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65f1f2ae84SBarry Smith }
66f1f2ae84SBarry Smith 
67d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPICommsDestroy(void)
68d71ae5a4SJacob Faibussowitsch {
69f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
70f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
71f1f2ae84SBarry Smith 
72f1f2ae84SBarry Smith   PetscFunctionBegin;
733ba16761SJacob Faibussowitsch   if (!PCMPICommSet) PetscFunctionReturn(PETSC_SUCCESS);
74f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
75f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
76f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
77f1f2ae84SBarry Smith     if (PCMPIComms[i] != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_free(&PCMPIComms[i]));
78f1f2ae84SBarry Smith   }
79f1f2ae84SBarry Smith   PCMPICommSet = PETSC_FALSE;
803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81f1f2ae84SBarry Smith }
82f1f2ae84SBarry Smith 
83d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPICreate(PC pc)
84d71ae5a4SJacob Faibussowitsch {
85f1f2ae84SBarry Smith   PC_MPI     *km   = pc ? (PC_MPI *)pc->data : NULL;
86f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
87f1f2ae84SBarry Smith   KSP         ksp;
88f1f2ae84SBarry Smith   PetscInt    N[2], mincntperrank = 0;
89f1f2ae84SBarry Smith   PetscMPIInt size;
90f1f2ae84SBarry Smith   Mat         sA;
913821be0aSBarry Smith   char       *cprefix = NULL;
92f1f2ae84SBarry Smith   PetscMPIInt len     = 0;
93f1f2ae84SBarry Smith 
94f1f2ae84SBarry Smith   PetscFunctionBegin;
95f1f2ae84SBarry Smith   if (!PCMPICommSet) PetscCall(PCMPICommsCreate());
96f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
97f1f2ae84SBarry Smith   if (pc) {
98f1f2ae84SBarry Smith     if (size == 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Warning: Running KSP type of MPI on a one rank MPI run, this will be less efficient then not using this type\n"));
99f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
100f1f2ae84SBarry Smith     PetscCall(MatGetSize(sA, &N[0], &N[1]));
101f1f2ae84SBarry Smith   }
102f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(N, 2, MPIU_INT, 0, comm));
103f1f2ae84SBarry Smith 
104f1f2ae84SBarry Smith   /* choose a suitable sized MPI_Comm for the problem to be solved on */
105f1f2ae84SBarry Smith   if (km) mincntperrank = km->mincntperrank;
106f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&mincntperrank, 1, MPI_INT, 0, comm));
107f1f2ae84SBarry Smith   comm = PCMPIComms[PetscMin(size, PetscMax(1, N[0] / mincntperrank)) - 1];
108f1f2ae84SBarry Smith   if (comm == MPI_COMM_NULL) {
109f1f2ae84SBarry Smith     ksp = NULL;
1103ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
111f1f2ae84SBarry Smith   }
1120316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
113f1f2ae84SBarry Smith   PetscCall(KSPCreate(comm, &ksp));
1143821be0aSBarry Smith   PetscCall(KSPSetNestLevel(ksp, 1));
1153821be0aSBarry Smith   PetscCall(PetscObjectSetTabLevel((PetscObject)ksp, 1));
1160316ec64SBarry Smith   PetscCall(PetscLogStagePop());
117f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gather(&ksp, 1, MPI_AINT, pc ? km->ksps : NULL, 1, MPI_AINT, 0, comm));
118f1f2ae84SBarry Smith   if (pc) {
119f1f2ae84SBarry Smith     size_t      slen;
1203821be0aSBarry Smith     const char *prefix = NULL;
121f9818f3cSJose E. Roman     char       *found  = NULL;
122f1f2ae84SBarry Smith 
123f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
124dad3da8eSBarry Smith     PCMPIKSPCounts[size - 1]++;
1253821be0aSBarry Smith     /* Created KSP gets prefix of PC minus the mpi_linear_solver_server_ portion */
1263821be0aSBarry Smith     PetscCall(PCGetOptionsPrefix(pc, &prefix));
1273821be0aSBarry Smith     PetscCheck(prefix, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing required prefix");
1283821be0aSBarry Smith     PetscCall(PetscStrallocpy(prefix, &cprefix));
1293821be0aSBarry Smith     PetscCall(PetscStrstr(cprefix, "mpi_linear_solver_server_", &found));
1303821be0aSBarry Smith     PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing mpi_linear_solver_server_ portion of prefix");
1313821be0aSBarry Smith     *found = 0;
1323821be0aSBarry Smith     PetscCall(PetscStrlen(cprefix, &slen));
133f1f2ae84SBarry Smith     len = (PetscMPIInt)slen;
134f1f2ae84SBarry Smith   }
135f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&len, 1, MPI_INT, 0, comm));
136f1f2ae84SBarry Smith   if (len) {
1373821be0aSBarry Smith     if (!pc) PetscCall(PetscMalloc1(len + 1, &cprefix));
1383821be0aSBarry Smith     PetscCallMPI(MPI_Bcast(cprefix, len + 1, MPI_CHAR, 0, comm));
1393821be0aSBarry Smith     PetscCall(KSPSetOptionsPrefix(ksp, cprefix));
140f1f2ae84SBarry Smith   }
1413821be0aSBarry Smith   PetscCall(PetscFree(cprefix));
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143f1f2ae84SBarry Smith }
144f1f2ae84SBarry Smith 
145d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISetMat(PC pc)
146d71ae5a4SJacob Faibussowitsch {
147f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
148f1f2ae84SBarry Smith   Mat                A;
14968a21331SBarry Smith   PetscInt           m, n, *ia, *ja, j, bs;
150f1f2ae84SBarry Smith   Mat                sA;
151f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
152f1f2ae84SBarry Smith   KSP                ksp;
153f1f2ae84SBarry Smith   PetscLayout        layout;
154f1f2ae84SBarry Smith   const PetscInt    *IA = NULL, *JA = NULL;
155f1f2ae84SBarry Smith   const PetscInt    *range;
156f1f2ae84SBarry Smith   PetscMPIInt       *NZ = NULL, sendcounti[PC_MPI_MAX_RANKS], displi[PC_MPI_MAX_RANKS], *NZdispl = NULL, nz, size, i;
157f1f2ae84SBarry Smith   PetscScalar       *a;
158f1f2ae84SBarry Smith   const PetscScalar *sa               = NULL;
1593821be0aSBarry Smith   PetscInt           matproperties[8] = {0, 0, 0, 0, 0, 0, 0, 0};
1603821be0aSBarry Smith   char              *cprefix;
161f1f2ae84SBarry Smith 
162f1f2ae84SBarry Smith   PetscFunctionBegin;
163f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
1643ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
165*45682376SBarry Smith   PetscCall(PetscLogEventBegin(EventServerDist, NULL, NULL, NULL, NULL));
166f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
167f1f2ae84SBarry Smith   if (pc) {
16868a21331SBarry Smith     PetscBool   isset, issymmetric, ishermitian, isspd, isstructurallysymmetric;
1693821be0aSBarry Smith     const char *prefix;
1703821be0aSBarry Smith     size_t      clen;
17168a21331SBarry Smith 
172f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
173dad3da8eSBarry Smith     PCMPIMatCounts[size - 1]++;
174f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
17568a21331SBarry Smith     PetscCall(MatGetSize(sA, &matproperties[0], &matproperties[1]));
176dd0d27b1SBarry Smith     PetscCall(MatGetBlockSize(sA, &bs));
17768a21331SBarry Smith     matproperties[2] = bs;
17868a21331SBarry Smith     PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric));
17968a21331SBarry Smith     matproperties[3] = !isset ? 0 : (issymmetric ? 1 : 2);
18068a21331SBarry Smith     PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian));
18168a21331SBarry Smith     matproperties[4] = !isset ? 0 : (ishermitian ? 1 : 2);
18268a21331SBarry Smith     PetscCall(MatIsSPDKnown(sA, &isset, &isspd));
18368a21331SBarry Smith     matproperties[5] = !isset ? 0 : (isspd ? 1 : 2);
18468a21331SBarry Smith     PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric));
18568a21331SBarry Smith     matproperties[6] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2);
1867a99bfcaSBarry Smith     /* Created Mat gets prefix of input Mat PLUS the mpi_linear_solver_server_ portion */
1873821be0aSBarry Smith     PetscCall(MatGetOptionsPrefix(sA, &prefix));
1883821be0aSBarry Smith     PetscCall(PetscStrallocpy(prefix, &cprefix));
1893821be0aSBarry Smith     PetscCall(PetscStrlen(cprefix, &clen));
1903821be0aSBarry Smith     matproperties[7] = (PetscInt)clen;
191f1f2ae84SBarry Smith   }
1923821be0aSBarry Smith   PetscCallMPI(MPI_Bcast(matproperties, PETSC_STATIC_ARRAY_LENGTH(matproperties), MPIU_INT, 0, comm));
193f1f2ae84SBarry Smith 
19468a21331SBarry Smith   /* determine ownership ranges of matrix columns */
195f1f2ae84SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
19668a21331SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2]));
19768a21331SBarry Smith   PetscCall(PetscLayoutSetSize(layout, matproperties[1]));
198f1f2ae84SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
199f1f2ae84SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &n));
20068a21331SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
20168a21331SBarry Smith 
20268a21331SBarry Smith   /* determine ownership ranges of matrix rows */
20368a21331SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
20468a21331SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2]));
20568a21331SBarry Smith   PetscCall(PetscLayoutSetSize(layout, matproperties[0]));
20668a21331SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
20768a21331SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &m));
208f1f2ae84SBarry Smith 
209f1f2ae84SBarry Smith   /* copy over the matrix nonzero structure and values */
210f1f2ae84SBarry Smith   if (pc) {
211f1f2ae84SBarry Smith     NZ      = km->NZ;
212f1f2ae84SBarry Smith     NZdispl = km->NZdispl;
213f1f2ae84SBarry Smith     PetscCall(PetscLayoutGetRanges(layout, &range));
214f1f2ae84SBarry Smith     PetscCall(MatGetRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
215f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
216f1f2ae84SBarry Smith       sendcounti[i] = (PetscMPIInt)(1 + range[i + 1] - range[i]);
217f1f2ae84SBarry Smith       NZ[i]         = (PetscMPIInt)(IA[range[i + 1]] - IA[range[i]]);
218f1f2ae84SBarry Smith     }
219f1f2ae84SBarry Smith     displi[0]  = 0;
220f1f2ae84SBarry Smith     NZdispl[0] = 0;
221f1f2ae84SBarry Smith     for (j = 1; j < size; j++) {
222f1f2ae84SBarry Smith       displi[j]  = displi[j - 1] + sendcounti[j - 1] - 1;
223f1f2ae84SBarry Smith       NZdispl[j] = NZdispl[j - 1] + NZ[j - 1];
224f1f2ae84SBarry Smith     }
225f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
226f1f2ae84SBarry Smith   }
227f1f2ae84SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
228f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(NZ, 1, MPI_INT, &nz, 1, MPI_INT, 0, comm));
229f1f2ae84SBarry Smith 
230f1f2ae84SBarry Smith   PetscCall(PetscMalloc3(n + 1, &ia, nz, &ja, nz, &a));
231f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(IA, sendcounti, displi, MPIU_INT, ia, n + 1, MPIU_INT, 0, comm));
232f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(JA, NZ, NZdispl, MPIU_INT, ja, nz, MPIU_INT, 0, comm));
233f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, NZ, NZdispl, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
234f1f2ae84SBarry Smith 
235f1f2ae84SBarry Smith   if (pc) {
236f1f2ae84SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
237f1f2ae84SBarry Smith     PetscCall(MatRestoreRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
238f1f2ae84SBarry Smith   }
239f1f2ae84SBarry Smith 
240ad540459SPierre Jolivet   for (j = 1; j < n + 1; j++) ia[j] -= ia[0];
241f1f2ae84SBarry Smith   ia[0] = 0;
2420316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
2433821be0aSBarry Smith   PetscCall(MatCreate(comm, &A));
2443821be0aSBarry Smith   if (matproperties[7] > 0) {
2453821be0aSBarry Smith     if (!pc) PetscCall(PetscMalloc1(matproperties[7] + 1, &cprefix));
2463821be0aSBarry Smith     PetscCallMPI(MPI_Bcast(cprefix, matproperties[7] + 1, MPI_CHAR, 0, comm));
2473821be0aSBarry Smith     PetscCall(MatSetOptionsPrefix(A, cprefix));
2483821be0aSBarry Smith     PetscCall(PetscFree(cprefix));
2493821be0aSBarry Smith   }
2503821be0aSBarry Smith   PetscCall(MatAppendOptionsPrefix(A, "mpi_linear_solver_server_"));
2513821be0aSBarry Smith   PetscCall(MatSetSizes(A, m, n, matproperties[0], matproperties[1]));
2523821be0aSBarry Smith   PetscCall(MatSetType(A, MATMPIAIJ));
2533821be0aSBarry Smith   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, a));
25468a21331SBarry Smith   PetscCall(MatSetBlockSize(A, matproperties[2]));
2553821be0aSBarry Smith 
25668a21331SBarry Smith   if (matproperties[3]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE));
25768a21331SBarry Smith   if (matproperties[4]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[4] == 1 ? PETSC_TRUE : PETSC_FALSE));
25868a21331SBarry Smith   if (matproperties[5]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[5] == 1 ? PETSC_TRUE : PETSC_FALSE));
25968a21331SBarry Smith   if (matproperties[6]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[6] == 1 ? PETSC_TRUE : PETSC_FALSE));
260f1f2ae84SBarry Smith 
261f1f2ae84SBarry Smith   PetscCall(PetscFree3(ia, ja, a));
262f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(ksp, A, A));
263f1f2ae84SBarry Smith   if (!ksp->vec_sol) PetscCall(MatCreateVecs(A, &ksp->vec_sol, &ksp->vec_rhs));
2640316ec64SBarry Smith   PetscCall(PetscLogStagePop());
265f1f2ae84SBarry Smith   if (pc) { /* needed for scatterv/gatherv of rhs and solution */
266f1f2ae84SBarry Smith     const PetscInt *range;
267f1f2ae84SBarry Smith 
268f1f2ae84SBarry Smith     PetscCall(VecGetOwnershipRanges(ksp->vec_sol, &range));
269f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
270f1f2ae84SBarry Smith       km->sendcount[i] = (PetscMPIInt)(range[i + 1] - range[i]);
271f1f2ae84SBarry Smith       km->displ[i]     = (PetscMPIInt)range[i];
272f1f2ae84SBarry Smith     }
273f1f2ae84SBarry Smith   }
274f1f2ae84SBarry Smith   PetscCall(MatDestroy(&A));
275*45682376SBarry Smith   PetscCall(PetscLogEventEnd(EventServerDist, NULL, NULL, NULL, NULL));
276f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(ksp));
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
278f1f2ae84SBarry Smith }
279f1f2ae84SBarry Smith 
280d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIUpdateMatValues(PC pc)
281d71ae5a4SJacob Faibussowitsch {
282f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
283f1f2ae84SBarry Smith   KSP                ksp;
284f1f2ae84SBarry Smith   Mat                sA, A;
285f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
286f1f2ae84SBarry Smith   PetscScalar       *a;
287f1f2ae84SBarry Smith   PetscCount         nz;
288f1f2ae84SBarry Smith   const PetscScalar *sa = NULL;
289dad3da8eSBarry Smith   PetscMPIInt        size;
29068a21331SBarry Smith   PetscInt           matproperties[4] = {0, 0, 0, 0};
291f1f2ae84SBarry Smith 
292f1f2ae84SBarry Smith   PetscFunctionBegin;
293f1f2ae84SBarry Smith   if (pc) {
294f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
295f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
296f1f2ae84SBarry Smith   }
297f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
2983ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
299*45682376SBarry Smith   PetscCall(PetscLogEventBegin(EventServerDist, NULL, NULL, NULL, NULL));
300f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
301dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
302dad3da8eSBarry Smith   PCMPIMatCounts[size - 1]++;
303f1f2ae84SBarry Smith   PetscCall(KSPGetOperators(ksp, NULL, &A));
304f1f2ae84SBarry Smith   PetscCall(MatMPIAIJGetNumberNonzeros(A, &nz));
305f1f2ae84SBarry Smith   PetscCall(PetscMalloc1(nz, &a));
306f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, pc ? km->NZ : NULL, pc ? km->NZdispl : NULL, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
30768a21331SBarry Smith   if (pc) {
30868a21331SBarry Smith     PetscBool isset, issymmetric, ishermitian, isspd, isstructurallysymmetric;
30968a21331SBarry Smith 
31068a21331SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
31168a21331SBarry Smith 
31268a21331SBarry Smith     PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric));
31368a21331SBarry Smith     matproperties[0] = !isset ? 0 : (issymmetric ? 1 : 2);
31468a21331SBarry Smith     PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian));
31568a21331SBarry Smith     matproperties[1] = !isset ? 0 : (ishermitian ? 1 : 2);
31668a21331SBarry Smith     PetscCall(MatIsSPDKnown(sA, &isset, &isspd));
31768a21331SBarry Smith     matproperties[2] = !isset ? 0 : (isspd ? 1 : 2);
31868a21331SBarry Smith     PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric));
31968a21331SBarry Smith     matproperties[3] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2);
32068a21331SBarry Smith   }
321f1f2ae84SBarry Smith   PetscCall(MatUpdateMPIAIJWithArray(A, a));
322f1f2ae84SBarry Smith   PetscCall(PetscFree(a));
32368a21331SBarry Smith   PetscCallMPI(MPI_Bcast(matproperties, 4, MPIU_INT, 0, comm));
32468a21331SBarry Smith   /* if any of these properties was previously set and is now not set this will result in incorrect properties in A since there is no way to unset a property */
32568a21331SBarry Smith   if (matproperties[0]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[0] == 1 ? PETSC_TRUE : PETSC_FALSE));
32668a21331SBarry Smith   if (matproperties[1]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[1] == 1 ? PETSC_TRUE : PETSC_FALSE));
32768a21331SBarry Smith   if (matproperties[2]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[2] == 1 ? PETSC_TRUE : PETSC_FALSE));
32868a21331SBarry Smith   if (matproperties[3]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE));
329*45682376SBarry Smith   PetscCall(PetscLogEventEnd(EventServerDist, NULL, NULL, NULL, NULL));
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331f1f2ae84SBarry Smith }
332f1f2ae84SBarry Smith 
333d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISolve(PC pc, Vec B, Vec X)
334d71ae5a4SJacob Faibussowitsch {
335f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
336f1f2ae84SBarry Smith   KSP                ksp;
337f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
338f1f2ae84SBarry Smith   const PetscScalar *sb   = NULL, *x;
339f1f2ae84SBarry Smith   PetscScalar       *b, *sx = NULL;
3405316cbedSBarry Smith   PetscInt           its, n;
3415316cbedSBarry Smith   PetscMPIInt        size;
342f1f2ae84SBarry Smith 
343f1f2ae84SBarry Smith   PetscFunctionBegin;
344f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : &ksp, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
3453ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
346*45682376SBarry Smith   PetscCall(PetscLogEventBegin(EventServerDist, NULL, NULL, NULL, NULL));
347f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
348f1f2ae84SBarry Smith 
349f1f2ae84SBarry Smith   /* TODO: optimize code to not require building counts/displ every time */
350f1f2ae84SBarry Smith 
351f1f2ae84SBarry Smith   /* scatterv rhs */
352dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
3535316cbedSBarry Smith   if (pc) {
3545316cbedSBarry Smith     PetscInt N;
3555316cbedSBarry Smith 
356dad3da8eSBarry Smith     PCMPISolveCounts[size - 1]++;
35768a21331SBarry Smith     PetscCall(MatGetSize(pc->pmat, &N, NULL));
3585316cbedSBarry Smith     PCMPISizes[size - 1] += N;
359f1f2ae84SBarry Smith     PetscCall(VecGetArrayRead(B, &sb));
360f1f2ae84SBarry Smith   }
361f1f2ae84SBarry Smith   PetscCall(VecGetLocalSize(ksp->vec_rhs, &n));
362f1f2ae84SBarry Smith   PetscCall(VecGetArray(ksp->vec_rhs, &b));
363f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sb, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, b, n, MPIU_SCALAR, 0, comm));
364f1f2ae84SBarry Smith   PetscCall(VecRestoreArray(ksp->vec_rhs, &b));
365f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArrayRead(B, &sb));
366f1f2ae84SBarry Smith 
367*45682376SBarry Smith   PetscCall(PetscLogEventEnd(EventServerDist, NULL, NULL, NULL, NULL));
3680316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
369f1f2ae84SBarry Smith   PetscCall(KSPSolve(ksp, NULL, NULL));
3700316ec64SBarry Smith   PetscCall(PetscLogStagePop());
371*45682376SBarry Smith   PetscCall(PetscLogEventBegin(EventServerDist, NULL, NULL, NULL, NULL));
3725316cbedSBarry Smith   PetscCall(KSPGetIterationNumber(ksp, &its));
3735316cbedSBarry Smith   PCMPIIterations[size - 1] += its;
374f1f2ae84SBarry Smith 
375f1f2ae84SBarry Smith   /* gather solution */
376f1f2ae84SBarry Smith   PetscCall(VecGetArrayRead(ksp->vec_sol, &x));
377f1f2ae84SBarry Smith   if (pc) PetscCall(VecGetArray(X, &sx));
378f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gatherv(x, n, MPIU_SCALAR, sx, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, 0, comm));
379f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArray(X, &sx));
380f1f2ae84SBarry Smith   PetscCall(VecRestoreArrayRead(ksp->vec_sol, &x));
381*45682376SBarry Smith   PetscCall(PetscLogEventEnd(EventServerDist, NULL, NULL, NULL, NULL));
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383f1f2ae84SBarry Smith }
384f1f2ae84SBarry Smith 
385d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIDestroy(PC pc)
386d71ae5a4SJacob Faibussowitsch {
387f1f2ae84SBarry Smith   PC_MPI  *km = pc ? (PC_MPI *)pc->data : NULL;
388f1f2ae84SBarry Smith   KSP      ksp;
389f1f2ae84SBarry Smith   MPI_Comm comm = PC_MPI_COMM_WORLD;
390f1f2ae84SBarry Smith 
391f1f2ae84SBarry Smith   PetscFunctionBegin;
392f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
3933ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
394c7d372c4SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
395f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&ksp));
396c7d372c4SBarry Smith   PetscCall(PetscLogStagePop());
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398f1f2ae84SBarry Smith }
399f1f2ae84SBarry Smith 
4003821be0aSBarry Smith PetscBool PCMPIServerActive = PETSC_FALSE;
4013821be0aSBarry Smith 
402f1f2ae84SBarry Smith /*@C
4037a99bfcaSBarry Smith   PCMPIServerBegin - starts a server that runs on the `rank != 0` MPI processes waiting to process requests for
404f1580f4eSBarry Smith   parallel `KSP` solves and management of parallel `KSP` objects.
405f1f2ae84SBarry Smith 
4063821be0aSBarry Smith   Logically Collective on all MPI processes except rank 0
407f1f2ae84SBarry Smith 
408f1580f4eSBarry Smith   Options Database Keys:
409f1f2ae84SBarry Smith + -mpi_linear_solver_server      - causes the PETSc program to start in MPI linear solver server mode where only the first MPI rank runs user code
4103821be0aSBarry Smith - -mpi_linear_solver_server_view - displays information about all the linear systems solved by the MPI linear solver server at the conclusion of the program
411f1f2ae84SBarry Smith 
41220f4b53cSBarry Smith   Level: developer
41320f4b53cSBarry Smith 
414f1580f4eSBarry Smith   Note:
415f1f2ae84SBarry Smith   This is normally started automatically in `PetscInitialize()` when the option is provided
416f1f2ae84SBarry Smith 
4173821be0aSBarry Smith   See `PCMPI` for information on using the solver with a `KSP` object
4183821be0aSBarry Smith 
419f1f2ae84SBarry Smith   Developer Notes:
4203821be0aSBarry Smith   When called on MPI rank 0 this sets `PETSC_COMM_WORLD` to `PETSC_COMM_SELF` to allow a main program
421f1580f4eSBarry Smith   written with `PETSC_COMM_WORLD` to run correctly on the single rank while all the ranks
422f1580f4eSBarry Smith   (that would normally be sharing `PETSC_COMM_WORLD`) to run the solver server.
423f1f2ae84SBarry Smith 
424f1580f4eSBarry Smith   Can this be integrated into the `PetscDevice` abstraction that is currently being developed?
425f1f2ae84SBarry Smith 
4263821be0aSBarry Smith   Conceivably `PCREDISTRIBUTE` could be organized in a similar manner to simplify its usage
4273821be0aSBarry Smith 
4283821be0aSBarry Smith   This could be implemented directly at the `KSP` level instead of using the `PCMPI` wrapper object
4293821be0aSBarry Smith 
430baca6076SPierre Jolivet   The code could be extended to allow an MPI + OpenMP application to use the linear solver server concept across all shared-memory
4313821be0aSBarry Smith   nodes with a single MPI process per node for the user application but multiple MPI processes per node for the linear solver.
4323821be0aSBarry Smith 
4337a99bfcaSBarry Smith   The concept could also be extended for users's callbacks for `SNES`, `TS`, and `Tao` where the `SNESSolve()` for example, runs on
4343821be0aSBarry Smith   all MPI processes but the user callback only runs on one MPI process per node.
4353821be0aSBarry Smith 
4363821be0aSBarry Smith   PETSc could also be extended with an MPI-less API that provides access to PETSc's solvers without any reference to MPI, essentially remove
4373821be0aSBarry Smith   the `MPI_Comm` argument from PETSc calls.
4383821be0aSBarry Smith 
4393821be0aSBarry Smith .seealso: [](sec_pcmpi), `PCMPIServerEnd()`, `PCMPI`, `KSPCheckPCMPI()`
440f1f2ae84SBarry Smith @*/
441d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerBegin(void)
442d71ae5a4SJacob Faibussowitsch {
443f1f2ae84SBarry Smith   PetscMPIInt rank;
444f1f2ae84SBarry Smith 
445f1f2ae84SBarry Smith   PetscFunctionBegin;
4469d3446b2SPierre Jolivet   PetscCall(PetscInfo(NULL, "Starting MPI Linear Solver Server\n"));
4475e1a0e3cSBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY)) {
4485e1a0e3cSBarry Smith     PetscCall(VecInitializePackage());
4495e1a0e3cSBarry Smith     PetscCall(MatInitializePackage());
4505e1a0e3cSBarry Smith     PetscCall(DMInitializePackage());
4515e1a0e3cSBarry Smith     PetscCall(PCInitializePackage());
4525e1a0e3cSBarry Smith     PetscCall(KSPInitializePackage());
4535e1a0e3cSBarry Smith     PetscCall(SNESInitializePackage());
4545e1a0e3cSBarry Smith     PetscCall(TSInitializePackage());
4555e1a0e3cSBarry Smith     PetscCall(TaoInitializePackage());
4565e1a0e3cSBarry Smith   }
457956255efSBarry Smith   PetscCall(PetscLogStageRegister("PCMPI", &PCMPIStage));
458*45682376SBarry Smith   PetscCall(PetscLogEventRegister("ServerDist", PC_CLASSID, &EventServerDist));
4595e1a0e3cSBarry Smith 
460f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(PC_MPI_COMM_WORLD, &rank));
461f1f2ae84SBarry Smith   if (rank == 0) {
462f1f2ae84SBarry Smith     PETSC_COMM_WORLD  = PETSC_COMM_SELF;
4633821be0aSBarry Smith     PCMPIServerActive = PETSC_TRUE;
4643ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
465f1f2ae84SBarry Smith   }
466f1f2ae84SBarry Smith 
467f1f2ae84SBarry Smith   while (PETSC_TRUE) {
468f1f2ae84SBarry Smith     PCMPICommand request = PCMPI_CREATE;
469f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
470f1f2ae84SBarry Smith     switch (request) {
471d71ae5a4SJacob Faibussowitsch     case PCMPI_CREATE:
472d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPICreate(NULL));
473d71ae5a4SJacob Faibussowitsch       break;
474d71ae5a4SJacob Faibussowitsch     case PCMPI_SET_MAT:
475d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPISetMat(NULL));
476d71ae5a4SJacob Faibussowitsch       break;
477d71ae5a4SJacob Faibussowitsch     case PCMPI_UPDATE_MAT_VALUES:
478d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPIUpdateMatValues(NULL));
479d71ae5a4SJacob Faibussowitsch       break;
480f1f2ae84SBarry Smith     case PCMPI_VIEW:
481f1f2ae84SBarry Smith       // PetscCall(PCMPIView(NULL));
482f1f2ae84SBarry Smith       break;
483d71ae5a4SJacob Faibussowitsch     case PCMPI_SOLVE:
484d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPISolve(NULL, NULL, NULL));
485d71ae5a4SJacob Faibussowitsch       break;
486d71ae5a4SJacob Faibussowitsch     case PCMPI_DESTROY:
487d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPIDestroy(NULL));
488d71ae5a4SJacob Faibussowitsch       break;
489f1f2ae84SBarry Smith     case PCMPI_EXIT:
490f1f2ae84SBarry Smith       PetscCall(PetscFinalize());
491f1f2ae84SBarry Smith       exit(0); /* not sure if this is a good idea, but cannot return because it will run users main program */
492f1f2ae84SBarry Smith       break;
493d71ae5a4SJacob Faibussowitsch     default:
494d71ae5a4SJacob Faibussowitsch       break;
495f1f2ae84SBarry Smith     }
496f1f2ae84SBarry Smith   }
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
498f1f2ae84SBarry Smith }
499f1f2ae84SBarry Smith 
500f1f2ae84SBarry Smith /*@C
501f1f2ae84SBarry Smith   PCMPIServerEnd - ends a server that runs on the rank != 0 MPI processes waiting to process requests for
502f1580f4eSBarry Smith   parallel KSP solves and management of parallel `KSP` objects.
503f1f2ae84SBarry Smith 
50420f4b53cSBarry Smith   Logically Collective on all MPI ranks except 0
50520f4b53cSBarry Smith 
50620f4b53cSBarry Smith   Level: developer
507f1f2ae84SBarry Smith 
508f1580f4eSBarry Smith   Note:
509f1f2ae84SBarry Smith   This is normally ended automatically in `PetscFinalize()` when the option is provided
510f1f2ae84SBarry Smith 
5113821be0aSBarry Smith .seealso: [](sec_pcmpi), `PCMPIServerBegin()`, `PCMPI`, `KSPCheckPCMPI()`
512f1f2ae84SBarry Smith @*/
513d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerEnd(void)
514d71ae5a4SJacob Faibussowitsch {
515f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_EXIT;
516f1f2ae84SBarry Smith 
517f1f2ae84SBarry Smith   PetscFunctionBegin;
518f1f2ae84SBarry Smith   if (PetscGlobalRank == 0) {
519f1f2ae84SBarry Smith     PetscViewer       viewer = NULL;
520f1f2ae84SBarry Smith     PetscViewerFormat format;
521f1f2ae84SBarry Smith 
522f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
523f1f2ae84SBarry Smith     PETSC_COMM_WORLD = MPI_COMM_WORLD; /* could use PC_MPI_COMM_WORLD */
524f1f2ae84SBarry Smith     PetscOptionsBegin(PETSC_COMM_SELF, NULL, "MPI linear solver server options", NULL);
525f1f2ae84SBarry Smith     PetscCall(PetscOptionsViewer("-mpi_linear_solver_server_view", "View information about system solved with the server", "PCMPI", &viewer, &format, NULL));
526f1f2ae84SBarry Smith     PetscOptionsEnd();
527f1f2ae84SBarry Smith     if (viewer) {
528f1f2ae84SBarry Smith       PetscBool isascii;
529f1f2ae84SBarry Smith 
530f1f2ae84SBarry Smith       PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
531f1f2ae84SBarry Smith       if (isascii) {
532f1f2ae84SBarry Smith         PetscMPIInt size;
5335316cbedSBarry Smith         PetscMPIInt i;
534f1f2ae84SBarry Smith 
535f1f2ae84SBarry Smith         PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
5365316cbedSBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "MPI linear solver server statistics:\n"));
5375316cbedSBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "    Ranks        KSPSolve()s     Mats        KSPs       Avg. Size      Avg. Its\n"));
5385316cbedSBarry Smith         if (PCMPIKSPCountsSeq) {
5395316cbedSBarry Smith           PetscCall(PetscViewerASCIIPrintf(viewer, "  Sequential         %" PetscInt_FMT "                         %" PetscInt_FMT "            %" PetscInt_FMT "           %" PetscInt_FMT "\n", PCMPISolveCountsSeq, PCMPIKSPCountsSeq, PCMPISizesSeq / PCMPISolveCountsSeq, PCMPIIterationsSeq / PCMPISolveCountsSeq));
540f1f2ae84SBarry Smith         }
5415316cbedSBarry Smith         for (i = 0; i < size; i++) {
5425316cbedSBarry Smith           if (PCMPIKSPCounts[i]) {
5435316cbedSBarry Smith             PetscCall(PetscViewerASCIIPrintf(viewer, "     %d               %" PetscInt_FMT "            %" PetscInt_FMT "           %" PetscInt_FMT "            %" PetscInt_FMT "            %" PetscInt_FMT "\n", i + 1, PCMPISolveCounts[i], PCMPIMatCounts[i], PCMPIKSPCounts[i], PCMPISizes[i] / PCMPISolveCounts[i], PCMPIIterations[i] / PCMPISolveCounts[i]));
5445316cbedSBarry Smith           }
5455316cbedSBarry Smith         }
546f1f2ae84SBarry Smith       }
547648c30bcSBarry Smith       PetscCall(PetscViewerDestroy(&viewer));
548f1f2ae84SBarry Smith     }
549f1f2ae84SBarry Smith   }
550f1f2ae84SBarry Smith   PetscCall(PCMPICommsDestroy());
5513821be0aSBarry Smith   PCMPIServerActive = PETSC_FALSE;
5523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
553f1f2ae84SBarry Smith }
554f1f2ae84SBarry Smith 
555f1f2ae84SBarry Smith /*
556f1f2ae84SBarry Smith     This version is used in the trivial case when the MPI parallel solver server is running on just the original MPI rank 0
557f1f2ae84SBarry Smith     because, for example, the problem is small. This version is more efficient because it does not require copying any data
558f1f2ae84SBarry Smith */
559d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Seq(PC pc)
560d71ae5a4SJacob Faibussowitsch {
561f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
562f1f2ae84SBarry Smith   Mat         sA;
563f1f2ae84SBarry Smith   const char *prefix;
564f9818f3cSJose E. Roman   char       *found = NULL, *cprefix;
565f1f2ae84SBarry Smith 
566f1f2ae84SBarry Smith   PetscFunctionBegin;
567f1f2ae84SBarry Smith   PetscCall(PCGetOperators(pc, NULL, &sA));
568f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
569f1f2ae84SBarry Smith   PetscCall(KSPCreate(PETSC_COMM_SELF, &km->ksps[0]));
5703821be0aSBarry Smith   PetscCall(KSPSetNestLevel(km->ksps[0], 1));
5713821be0aSBarry Smith   PetscCall(PetscObjectSetTabLevel((PetscObject)km->ksps[0], 1));
5723821be0aSBarry Smith 
5733821be0aSBarry Smith   /* Created KSP gets prefix of PC minus the mpi_linear_solver_server_ portion */
5743821be0aSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
5753821be0aSBarry Smith   PetscCheck(prefix, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing required prefix");
5763821be0aSBarry Smith   PetscCall(PetscStrallocpy(prefix, &cprefix));
5773821be0aSBarry Smith   PetscCall(PetscStrstr(cprefix, "mpi_linear_solver_server_", &found));
5783821be0aSBarry Smith   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing mpi_linear_solver_server_ portion of prefix");
5793821be0aSBarry Smith   *found = 0;
5803821be0aSBarry Smith   PetscCall(KSPSetOptionsPrefix(km->ksps[0], cprefix));
5813821be0aSBarry Smith   PetscCall(PetscFree(cprefix));
5823821be0aSBarry Smith 
583f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(km->ksps[0], sA, sA));
584f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(km->ksps[0]));
585f1f2ae84SBarry Smith   PetscCall(KSPSetUp(km->ksps[0]));
5863ba16761SJacob Faibussowitsch   PetscCall(PetscInfo((PetscObject)pc, "MPI parallel linear solver system is being solved directly on rank 0 due to its small size\n"));
587f1f2ae84SBarry Smith   PCMPIKSPCountsSeq++;
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
589f1f2ae84SBarry Smith }
590f1f2ae84SBarry Smith 
591d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Seq(PC pc, Vec b, Vec x)
592d71ae5a4SJacob Faibussowitsch {
593f1f2ae84SBarry Smith   PC_MPI  *km = (PC_MPI *)pc->data;
5945316cbedSBarry Smith   PetscInt its, n;
5955316cbedSBarry Smith   Mat      A;
596f1f2ae84SBarry Smith 
597f1f2ae84SBarry Smith   PetscFunctionBegin;
598f1f2ae84SBarry Smith   PetscCall(KSPSolve(km->ksps[0], b, x));
5995316cbedSBarry Smith   PetscCall(KSPGetIterationNumber(km->ksps[0], &its));
600f1f2ae84SBarry Smith   PCMPISolveCountsSeq++;
6015316cbedSBarry Smith   PCMPIIterationsSeq += its;
6025316cbedSBarry Smith   PetscCall(KSPGetOperators(km->ksps[0], NULL, &A));
6035316cbedSBarry Smith   PetscCall(MatGetSize(A, &n, NULL));
6045316cbedSBarry Smith   PCMPISizesSeq += n;
6053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
606f1f2ae84SBarry Smith }
607f1f2ae84SBarry Smith 
608d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Seq(PC pc, PetscViewer viewer)
609d71ae5a4SJacob Faibussowitsch {
610f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
611f1f2ae84SBarry Smith 
612f1f2ae84SBarry Smith   PetscFunctionBegin;
613f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Running MPI linear solver server directly on rank 0 due to its small size\n"));
614f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
6155316cbedSBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n"));
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617f1f2ae84SBarry Smith }
618f1f2ae84SBarry Smith 
619d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Seq(PC pc)
620d71ae5a4SJacob Faibussowitsch {
621f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
622f1f2ae84SBarry Smith 
623f1f2ae84SBarry Smith   PetscFunctionBegin;
624f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&km->ksps[0]));
625f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
6263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
627f1f2ae84SBarry Smith }
628f1f2ae84SBarry Smith 
629f1f2ae84SBarry Smith /*
630f1f2ae84SBarry Smith      PCSetUp_MPI - Trigger the creation of the MPI parallel PC and copy parts of the matrix and
631dd8e379bSPierre Jolivet      right-hand side to the parallel PC
632f1f2ae84SBarry Smith */
633d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_MPI(PC pc)
634d71ae5a4SJacob Faibussowitsch {
635f1f2ae84SBarry Smith   PC_MPI      *km = (PC_MPI *)pc->data;
636f1f2ae84SBarry Smith   PetscMPIInt  rank, size;
637f1f2ae84SBarry Smith   PCMPICommand request;
638f1f2ae84SBarry Smith   PetscBool    newmatrix = PETSC_FALSE;
639f1f2ae84SBarry Smith 
640f1f2ae84SBarry Smith   PetscFunctionBegin;
641f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
642f1f2ae84SBarry Smith   PetscCheck(rank == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PCMPI can only be used from 0th rank of MPI_COMM_WORLD. Perhaps a missing -mpi_linear_solver_server?");
643f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
644f1f2ae84SBarry Smith 
645f1f2ae84SBarry Smith   if (!pc->setupcalled) {
646f1f2ae84SBarry Smith     if (!km->alwaysuseserver) {
647f1f2ae84SBarry Smith       PetscInt n;
648f1f2ae84SBarry Smith       Mat      sA;
649f1f2ae84SBarry Smith       /* short circuit for small systems */
650f1f2ae84SBarry Smith       PetscCall(PCGetOperators(pc, &sA, &sA));
651f1f2ae84SBarry Smith       PetscCall(MatGetSize(sA, &n, NULL));
652f1f2ae84SBarry Smith       if (n < 2 * km->mincntperrank - 1 || size == 1) {
653f1f2ae84SBarry Smith         pc->ops->setup   = NULL;
654f1f2ae84SBarry Smith         pc->ops->apply   = PCApply_Seq;
655f1f2ae84SBarry Smith         pc->ops->destroy = PCDestroy_Seq;
656f1f2ae84SBarry Smith         pc->ops->view    = PCView_Seq;
657f1f2ae84SBarry Smith         PetscCall(PCSetUp_Seq(pc));
6583ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
659f1f2ae84SBarry Smith       }
660f1f2ae84SBarry Smith     }
661f1f2ae84SBarry Smith 
662f1f2ae84SBarry Smith     request = PCMPI_CREATE;
663f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
664f1f2ae84SBarry Smith     PetscCall(PCMPICreate(pc));
665f1f2ae84SBarry Smith     newmatrix = PETSC_TRUE;
6669371c9d4SSatish Balay   }
6679371c9d4SSatish Balay   if (pc->flag == DIFFERENT_NONZERO_PATTERN) newmatrix = PETSC_TRUE;
668f1f2ae84SBarry Smith 
669f1f2ae84SBarry Smith   if (newmatrix) {
6703ba16761SJacob Faibussowitsch     PetscCall(PetscInfo((PetscObject)pc, "New matrix or matrix has changed nonzero structure\n"));
671f1f2ae84SBarry Smith     request = PCMPI_SET_MAT;
672f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
673f1f2ae84SBarry Smith     PetscCall(PCMPISetMat(pc));
674f1f2ae84SBarry Smith   } else {
675bbea24aaSStefano Zampini     PetscCall(PetscInfo((PetscObject)pc, "Matrix has only changed nonzero values\n"));
676f1f2ae84SBarry Smith     request = PCMPI_UPDATE_MAT_VALUES;
677f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
678f1f2ae84SBarry Smith     PetscCall(PCMPIUpdateMatValues(pc));
679f1f2ae84SBarry Smith   }
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
681f1f2ae84SBarry Smith }
682f1f2ae84SBarry Smith 
683d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MPI(PC pc, Vec b, Vec x)
684d71ae5a4SJacob Faibussowitsch {
685f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_SOLVE;
686f1f2ae84SBarry Smith 
687f1f2ae84SBarry Smith   PetscFunctionBegin;
688f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
689f1f2ae84SBarry Smith   PetscCall(PCMPISolve(pc, b, x));
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
691f1f2ae84SBarry Smith }
692f1f2ae84SBarry Smith 
69366976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_MPI(PC pc)
694d71ae5a4SJacob Faibussowitsch {
695f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_DESTROY;
696f1f2ae84SBarry Smith 
697f1f2ae84SBarry Smith   PetscFunctionBegin;
698f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
699f1f2ae84SBarry Smith   PetscCall(PCMPIDestroy(pc));
700f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
7013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
702f1f2ae84SBarry Smith }
703f1f2ae84SBarry Smith 
704f1f2ae84SBarry Smith /*
705f1f2ae84SBarry Smith      PCView_MPI - Cannot call view on the MPI parallel KSP because other ranks do not have access to the viewer
706f1f2ae84SBarry Smith */
70766976f2fSJacob Faibussowitsch static PetscErrorCode PCView_MPI(PC pc, PetscViewer viewer)
708d71ae5a4SJacob Faibussowitsch {
709f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
710f1f2ae84SBarry Smith   MPI_Comm    comm;
711f1f2ae84SBarry Smith   PetscMPIInt size;
712f1f2ae84SBarry Smith 
713f1f2ae84SBarry Smith   PetscFunctionBegin;
714f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)km->ksps[0], &comm));
715f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
716f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Size of MPI communicator used for MPI parallel KSP solve %d\n", size));
717f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
71868a21331SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n"));
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
720f1f2ae84SBarry Smith }
721f1f2ae84SBarry Smith 
72266976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_MPI(PC pc, PetscOptionItems *PetscOptionsObject)
723d71ae5a4SJacob Faibussowitsch {
724f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
725f1f2ae84SBarry Smith 
726f1f2ae84SBarry Smith   PetscFunctionBegin;
727f1f2ae84SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPI linear solver server options");
7283821be0aSBarry Smith   PetscCall(PetscOptionsInt("-minimum_count_per_rank", "Desired minimum number of nonzeros per rank", "None", km->mincntperrank, &km->mincntperrank, NULL));
7293821be0aSBarry Smith   PetscCall(PetscOptionsBool("-always_use_server", "Use the server even if only one rank is used for the solve (for debugging)", "None", km->alwaysuseserver, &km->alwaysuseserver, NULL));
730f1f2ae84SBarry Smith   PetscOptionsHeadEnd();
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
732f1f2ae84SBarry Smith }
733f1f2ae84SBarry Smith 
734f1f2ae84SBarry Smith /*MC
735f1580f4eSBarry Smith      PCMPI - Calls an MPI parallel `KSP` to solve a linear system from user code running on one process
736f1f2ae84SBarry Smith 
7373821be0aSBarry Smith    Options Database Keys for the Server:
738f1f2ae84SBarry Smith +  -mpi_linear_solver_server - causes the PETSc program to start in MPI linear solver server mode where only the first MPI rank runs user code
7393821be0aSBarry Smith -  -mpi_linear_solver_server_view - displays information about all the linear systems solved by the MPI linear solver server
740f1f2ae84SBarry Smith 
7413821be0aSBarry Smith    Options Database Keys for a specific `KSP` object
7423821be0aSBarry Smith +  -[any_ksp_prefix]_mpi_linear_solver_server_minimum_count_per_rank - sets the minimum size of the linear system per MPI rank that the solver will strive for
7433821be0aSBarry Smith -  -[any_ksp_prefix]_mpi_linear_solver_server_always_use_server - use the server solver code even if the particular system is only solved on the process (for debugging and testing purposes)
7443821be0aSBarry Smith 
7453821be0aSBarry Smith    Level: developer
74620f4b53cSBarry Smith 
747f1f2ae84SBarry Smith    Notes:
74846bbbc36SPierre Jolivet    The options database prefix for the actual solver is any prefix provided before use to the original `KSP` with `KSPSetOptionsPrefix()`, mostly commonly no prefix is used.
749f1f2ae84SBarry Smith 
750f1f2ae84SBarry Smith    It can be particularly useful for user OpenMP code or potentially user GPU code.
751f1f2ae84SBarry Smith 
752dd8e379bSPierre Jolivet    When the program is running with a single MPI process then it directly uses the provided matrix and right-hand side
7533821be0aSBarry Smith    and does not need to distribute the matrix and vector to the various MPI processes; thus it incurs no extra overhead over just using the `KSP` directly.
754f1f2ae84SBarry Smith 
7553821be0aSBarry Smith    The solver options for actual solving `KSP` and `PC` must be controlled via the options database, calls to set options directly on the user level `KSP` and `PC` have no effect
7560316ec64SBarry Smith    because they are not the actual solver objects.
7570316ec64SBarry Smith 
7580316ec64SBarry Smith    When `-log_view` is used with this solver the events within the parallel solve are logging in their own stage. Some of the logging in the other
75968a21331SBarry Smith    stages will be confusing since the event times are only recorded on the 0th MPI rank, thus the percent of time in the events will be misleading.
7600316ec64SBarry Smith 
7613821be0aSBarry Smith    Developer Note:
7623821be0aSBarry Smith    This `PCType` is never directly selected by the user, it is set when the option `-mpi_linear_solver_server` is used and the `PC` is at the outer most nesting of
7633821be0aSBarry Smith    a `KSP`. The outer most `KSP` object is automatically set to `KSPPREONLY` and thus is not directly visible to the user.
7643821be0aSBarry Smith 
7653821be0aSBarry Smith .seealso: [](sec_pcmpi), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `PC`, `PCMPIServerBegin()`, `PCMPIServerEnd()`, `KSPCheckPCMPI()`
766f1f2ae84SBarry Smith M*/
767d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MPI(PC pc)
768d71ae5a4SJacob Faibussowitsch {
769f1f2ae84SBarry Smith   PC_MPI *km;
770f9818f3cSJose E. Roman   char   *found = NULL;
771f1f2ae84SBarry Smith 
772f1f2ae84SBarry Smith   PetscFunctionBegin;
7733821be0aSBarry Smith   PetscCall(PetscStrstr(((PetscObject)pc)->prefix, "mpi_linear_solver_server_", &found));
7743821be0aSBarry Smith   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI object prefix does not have mpi_linear_solver_server_");
7753821be0aSBarry Smith 
7763821be0aSBarry Smith   /* material from PCSetType() */
7773821be0aSBarry Smith   PetscTryTypeMethod(pc, destroy);
7783821be0aSBarry Smith   pc->ops->destroy = NULL;
7793821be0aSBarry Smith   pc->data         = NULL;
7803821be0aSBarry Smith 
7813821be0aSBarry Smith   PetscCall(PetscFunctionListDestroy(&((PetscObject)pc)->qlist));
7823821be0aSBarry Smith   PetscCall(PetscMemzero(pc->ops, sizeof(struct _PCOps)));
7833821be0aSBarry Smith   pc->modifysubmatrices  = NULL;
7843821be0aSBarry Smith   pc->modifysubmatricesP = NULL;
7853821be0aSBarry Smith   pc->setupcalled        = 0;
7863821be0aSBarry Smith 
7874dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&km));
788f1f2ae84SBarry Smith   pc->data = (void *)km;
789f1f2ae84SBarry Smith 
790f1f2ae84SBarry Smith   km->mincntperrank = 10000;
791f1f2ae84SBarry Smith 
792f1f2ae84SBarry Smith   pc->ops->setup          = PCSetUp_MPI;
793f1f2ae84SBarry Smith   pc->ops->apply          = PCApply_MPI;
794f1f2ae84SBarry Smith   pc->ops->destroy        = PCDestroy_MPI;
795f1f2ae84SBarry Smith   pc->ops->view           = PCView_MPI;
796f1f2ae84SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MPI;
7973821be0aSBarry Smith   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCMPI));
7983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
799f1f2ae84SBarry Smith }
800