xref: /petsc/src/ksp/pc/impls/mpi/pcmpi.c (revision 789afff49e22c0b31a41d280c30e86915c8d497b)
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>
15*789afff4SPierre Jolivet #include <petscts.h>
16*789afff4SPierre 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;
43f1f2ae84SBarry Smith 
44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPICommsCreate(void)
45d71ae5a4SJacob Faibussowitsch {
46f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
47f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
48f1f2ae84SBarry Smith 
49f1f2ae84SBarry Smith   PetscFunctionBegin;
50f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
51f1f2ae84SBarry 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");
52f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
53f1f2ae84SBarry Smith   /* comm for size 1 is useful only for debugging */
54f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
55f1f2ae84SBarry Smith     PetscMPIInt color = rank < i + 1 ? 0 : MPI_UNDEFINED;
56f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_split(comm, color, 0, &PCMPIComms[i]));
57f1f2ae84SBarry Smith     PCMPISolveCounts[i] = 0;
58f1f2ae84SBarry Smith     PCMPIKSPCounts[i]   = 0;
595316cbedSBarry Smith     PCMPIIterations[i]  = 0;
605316cbedSBarry Smith     PCMPISizes[i]       = 0;
61f1f2ae84SBarry Smith   }
62f1f2ae84SBarry Smith   PCMPICommSet = PETSC_TRUE;
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64f1f2ae84SBarry Smith }
65f1f2ae84SBarry Smith 
66d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPICommsDestroy(void)
67d71ae5a4SJacob Faibussowitsch {
68f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
69f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
70f1f2ae84SBarry Smith 
71f1f2ae84SBarry Smith   PetscFunctionBegin;
723ba16761SJacob Faibussowitsch   if (!PCMPICommSet) PetscFunctionReturn(PETSC_SUCCESS);
73f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
74f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
75f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
76f1f2ae84SBarry Smith     if (PCMPIComms[i] != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_free(&PCMPIComms[i]));
77f1f2ae84SBarry Smith   }
78f1f2ae84SBarry Smith   PCMPICommSet = PETSC_FALSE;
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80f1f2ae84SBarry Smith }
81f1f2ae84SBarry Smith 
82d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPICreate(PC pc)
83d71ae5a4SJacob Faibussowitsch {
84f1f2ae84SBarry Smith   PC_MPI     *km   = pc ? (PC_MPI *)pc->data : NULL;
85f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
86f1f2ae84SBarry Smith   KSP         ksp;
87f1f2ae84SBarry Smith   PetscInt    N[2], mincntperrank = 0;
88f1f2ae84SBarry Smith   PetscMPIInt size;
89f1f2ae84SBarry Smith   Mat         sA;
903821be0aSBarry Smith   char       *cprefix = NULL;
91f1f2ae84SBarry Smith   PetscMPIInt len     = 0;
92f1f2ae84SBarry Smith 
93f1f2ae84SBarry Smith   PetscFunctionBegin;
94f1f2ae84SBarry Smith   if (!PCMPICommSet) PetscCall(PCMPICommsCreate());
95f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
96f1f2ae84SBarry Smith   if (pc) {
97f1f2ae84SBarry 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"));
98f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
99f1f2ae84SBarry Smith     PetscCall(MatGetSize(sA, &N[0], &N[1]));
100f1f2ae84SBarry Smith   }
101f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(N, 2, MPIU_INT, 0, comm));
102f1f2ae84SBarry Smith 
103f1f2ae84SBarry Smith   /* choose a suitable sized MPI_Comm for the problem to be solved on */
104f1f2ae84SBarry Smith   if (km) mincntperrank = km->mincntperrank;
105f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&mincntperrank, 1, MPI_INT, 0, comm));
106f1f2ae84SBarry Smith   comm = PCMPIComms[PetscMin(size, PetscMax(1, N[0] / mincntperrank)) - 1];
107f1f2ae84SBarry Smith   if (comm == MPI_COMM_NULL) {
108f1f2ae84SBarry Smith     ksp = NULL;
1093ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
110f1f2ae84SBarry Smith   }
1110316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
112f1f2ae84SBarry Smith   PetscCall(KSPCreate(comm, &ksp));
1133821be0aSBarry Smith   PetscCall(KSPSetNestLevel(ksp, 1));
1143821be0aSBarry Smith   PetscCall(PetscObjectSetTabLevel((PetscObject)ksp, 1));
1150316ec64SBarry Smith   PetscCall(PetscLogStagePop());
116f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gather(&ksp, 1, MPI_AINT, pc ? km->ksps : NULL, 1, MPI_AINT, 0, comm));
117f1f2ae84SBarry Smith   if (pc) {
118f1f2ae84SBarry Smith     size_t      slen;
1193821be0aSBarry Smith     const char *prefix = NULL;
120f9818f3cSJose E. Roman     char       *found  = NULL;
121f1f2ae84SBarry Smith 
122f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
123dad3da8eSBarry Smith     PCMPIKSPCounts[size - 1]++;
1243821be0aSBarry Smith     /* Created KSP gets prefix of PC minus the mpi_linear_solver_server_ portion */
1253821be0aSBarry Smith     PetscCall(PCGetOptionsPrefix(pc, &prefix));
1263821be0aSBarry Smith     PetscCheck(prefix, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing required prefix");
1273821be0aSBarry Smith     PetscCall(PetscStrallocpy(prefix, &cprefix));
1283821be0aSBarry Smith     PetscCall(PetscStrstr(cprefix, "mpi_linear_solver_server_", &found));
1293821be0aSBarry Smith     PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing mpi_linear_solver_server_ portion of prefix");
1303821be0aSBarry Smith     *found = 0;
1313821be0aSBarry Smith     PetscCall(PetscStrlen(cprefix, &slen));
132f1f2ae84SBarry Smith     len = (PetscMPIInt)slen;
133f1f2ae84SBarry Smith   }
134f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&len, 1, MPI_INT, 0, comm));
135f1f2ae84SBarry Smith   if (len) {
1363821be0aSBarry Smith     if (!pc) PetscCall(PetscMalloc1(len + 1, &cprefix));
1373821be0aSBarry Smith     PetscCallMPI(MPI_Bcast(cprefix, len + 1, MPI_CHAR, 0, comm));
1383821be0aSBarry Smith     PetscCall(KSPSetOptionsPrefix(ksp, cprefix));
139f1f2ae84SBarry Smith   }
1403821be0aSBarry Smith   PetscCall(PetscFree(cprefix));
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142f1f2ae84SBarry Smith }
143f1f2ae84SBarry Smith 
144d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISetMat(PC pc)
145d71ae5a4SJacob Faibussowitsch {
146f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
147f1f2ae84SBarry Smith   Mat                A;
14868a21331SBarry Smith   PetscInt           m, n, *ia, *ja, j, bs;
149f1f2ae84SBarry Smith   Mat                sA;
150f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
151f1f2ae84SBarry Smith   KSP                ksp;
152f1f2ae84SBarry Smith   PetscLayout        layout;
153f1f2ae84SBarry Smith   const PetscInt    *IA = NULL, *JA = NULL;
154f1f2ae84SBarry Smith   const PetscInt    *range;
155f1f2ae84SBarry Smith   PetscMPIInt       *NZ = NULL, sendcounti[PC_MPI_MAX_RANKS], displi[PC_MPI_MAX_RANKS], *NZdispl = NULL, nz, size, i;
156f1f2ae84SBarry Smith   PetscScalar       *a;
157f1f2ae84SBarry Smith   const PetscScalar *sa               = NULL;
1583821be0aSBarry Smith   PetscInt           matproperties[8] = {0, 0, 0, 0, 0, 0, 0, 0};
1593821be0aSBarry Smith   char              *cprefix;
160f1f2ae84SBarry Smith 
161f1f2ae84SBarry Smith   PetscFunctionBegin;
162f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
1633ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
164f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
165f1f2ae84SBarry Smith   if (pc) {
16668a21331SBarry Smith     PetscBool   isset, issymmetric, ishermitian, isspd, isstructurallysymmetric;
1673821be0aSBarry Smith     const char *prefix;
1683821be0aSBarry Smith     size_t      clen;
16968a21331SBarry Smith 
170f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
171dad3da8eSBarry Smith     PCMPIMatCounts[size - 1]++;
172f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
17368a21331SBarry Smith     PetscCall(MatGetSize(sA, &matproperties[0], &matproperties[1]));
174dd0d27b1SBarry Smith     PetscCall(MatGetBlockSize(sA, &bs));
17568a21331SBarry Smith     matproperties[2] = bs;
17668a21331SBarry Smith     PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric));
17768a21331SBarry Smith     matproperties[3] = !isset ? 0 : (issymmetric ? 1 : 2);
17868a21331SBarry Smith     PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian));
17968a21331SBarry Smith     matproperties[4] = !isset ? 0 : (ishermitian ? 1 : 2);
18068a21331SBarry Smith     PetscCall(MatIsSPDKnown(sA, &isset, &isspd));
18168a21331SBarry Smith     matproperties[5] = !isset ? 0 : (isspd ? 1 : 2);
18268a21331SBarry Smith     PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric));
18368a21331SBarry Smith     matproperties[6] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2);
1847a99bfcaSBarry Smith     /* Created Mat gets prefix of input Mat PLUS the mpi_linear_solver_server_ portion */
1853821be0aSBarry Smith     PetscCall(MatGetOptionsPrefix(sA, &prefix));
1863821be0aSBarry Smith     PetscCall(PetscStrallocpy(prefix, &cprefix));
1873821be0aSBarry Smith     PetscCall(PetscStrlen(cprefix, &clen));
1883821be0aSBarry Smith     matproperties[7] = (PetscInt)clen;
189f1f2ae84SBarry Smith   }
1903821be0aSBarry Smith   PetscCallMPI(MPI_Bcast(matproperties, PETSC_STATIC_ARRAY_LENGTH(matproperties), MPIU_INT, 0, comm));
191f1f2ae84SBarry Smith 
19268a21331SBarry Smith   /* determine ownership ranges of matrix columns */
193f1f2ae84SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
19468a21331SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2]));
19568a21331SBarry Smith   PetscCall(PetscLayoutSetSize(layout, matproperties[1]));
196f1f2ae84SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
197f1f2ae84SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &n));
19868a21331SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
19968a21331SBarry Smith 
20068a21331SBarry Smith   /* determine ownership ranges of matrix rows */
20168a21331SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
20268a21331SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2]));
20368a21331SBarry Smith   PetscCall(PetscLayoutSetSize(layout, matproperties[0]));
20468a21331SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
20568a21331SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &m));
206f1f2ae84SBarry Smith 
207f1f2ae84SBarry Smith   /* copy over the matrix nonzero structure and values */
208f1f2ae84SBarry Smith   if (pc) {
209f1f2ae84SBarry Smith     NZ      = km->NZ;
210f1f2ae84SBarry Smith     NZdispl = km->NZdispl;
211f1f2ae84SBarry Smith     PetscCall(PetscLayoutGetRanges(layout, &range));
212f1f2ae84SBarry Smith     PetscCall(MatGetRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
213f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
214f1f2ae84SBarry Smith       sendcounti[i] = (PetscMPIInt)(1 + range[i + 1] - range[i]);
215f1f2ae84SBarry Smith       NZ[i]         = (PetscMPIInt)(IA[range[i + 1]] - IA[range[i]]);
216f1f2ae84SBarry Smith     }
217f1f2ae84SBarry Smith     displi[0]  = 0;
218f1f2ae84SBarry Smith     NZdispl[0] = 0;
219f1f2ae84SBarry Smith     for (j = 1; j < size; j++) {
220f1f2ae84SBarry Smith       displi[j]  = displi[j - 1] + sendcounti[j - 1] - 1;
221f1f2ae84SBarry Smith       NZdispl[j] = NZdispl[j - 1] + NZ[j - 1];
222f1f2ae84SBarry Smith     }
223f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
224f1f2ae84SBarry Smith   }
225f1f2ae84SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
226f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(NZ, 1, MPI_INT, &nz, 1, MPI_INT, 0, comm));
227f1f2ae84SBarry Smith 
228f1f2ae84SBarry Smith   PetscCall(PetscMalloc3(n + 1, &ia, nz, &ja, nz, &a));
229f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(IA, sendcounti, displi, MPIU_INT, ia, n + 1, MPIU_INT, 0, comm));
230f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(JA, NZ, NZdispl, MPIU_INT, ja, nz, MPIU_INT, 0, comm));
231f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, NZ, NZdispl, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
232f1f2ae84SBarry Smith 
233f1f2ae84SBarry Smith   if (pc) {
234f1f2ae84SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
235f1f2ae84SBarry Smith     PetscCall(MatRestoreRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
236f1f2ae84SBarry Smith   }
237f1f2ae84SBarry Smith 
238ad540459SPierre Jolivet   for (j = 1; j < n + 1; j++) ia[j] -= ia[0];
239f1f2ae84SBarry Smith   ia[0] = 0;
2400316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
2413821be0aSBarry Smith   PetscCall(MatCreate(comm, &A));
2423821be0aSBarry Smith   if (matproperties[7] > 0) {
2433821be0aSBarry Smith     if (!pc) PetscCall(PetscMalloc1(matproperties[7] + 1, &cprefix));
2443821be0aSBarry Smith     PetscCallMPI(MPI_Bcast(cprefix, matproperties[7] + 1, MPI_CHAR, 0, comm));
2453821be0aSBarry Smith     PetscCall(MatSetOptionsPrefix(A, cprefix));
2463821be0aSBarry Smith     PetscCall(PetscFree(cprefix));
2473821be0aSBarry Smith   }
2483821be0aSBarry Smith   PetscCall(MatAppendOptionsPrefix(A, "mpi_linear_solver_server_"));
2493821be0aSBarry Smith   PetscCall(MatSetSizes(A, m, n, matproperties[0], matproperties[1]));
2503821be0aSBarry Smith   PetscCall(MatSetType(A, MATMPIAIJ));
2513821be0aSBarry Smith   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, a));
25268a21331SBarry Smith   PetscCall(MatSetBlockSize(A, matproperties[2]));
2533821be0aSBarry Smith 
25468a21331SBarry Smith   if (matproperties[3]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE));
25568a21331SBarry Smith   if (matproperties[4]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[4] == 1 ? PETSC_TRUE : PETSC_FALSE));
25668a21331SBarry Smith   if (matproperties[5]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[5] == 1 ? PETSC_TRUE : PETSC_FALSE));
25768a21331SBarry Smith   if (matproperties[6]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[6] == 1 ? PETSC_TRUE : PETSC_FALSE));
258f1f2ae84SBarry Smith 
259f1f2ae84SBarry Smith   PetscCall(PetscFree3(ia, ja, a));
260f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(ksp, A, A));
261f1f2ae84SBarry Smith   if (!ksp->vec_sol) PetscCall(MatCreateVecs(A, &ksp->vec_sol, &ksp->vec_rhs));
2620316ec64SBarry Smith   PetscCall(PetscLogStagePop());
263f1f2ae84SBarry Smith   if (pc) { /* needed for scatterv/gatherv of rhs and solution */
264f1f2ae84SBarry Smith     const PetscInt *range;
265f1f2ae84SBarry Smith 
266f1f2ae84SBarry Smith     PetscCall(VecGetOwnershipRanges(ksp->vec_sol, &range));
267f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
268f1f2ae84SBarry Smith       km->sendcount[i] = (PetscMPIInt)(range[i + 1] - range[i]);
269f1f2ae84SBarry Smith       km->displ[i]     = (PetscMPIInt)range[i];
270f1f2ae84SBarry Smith     }
271f1f2ae84SBarry Smith   }
272f1f2ae84SBarry Smith   PetscCall(MatDestroy(&A));
273f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(ksp));
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275f1f2ae84SBarry Smith }
276f1f2ae84SBarry Smith 
277d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIUpdateMatValues(PC pc)
278d71ae5a4SJacob Faibussowitsch {
279f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
280f1f2ae84SBarry Smith   KSP                ksp;
281f1f2ae84SBarry Smith   Mat                sA, A;
282f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
283f1f2ae84SBarry Smith   PetscScalar       *a;
284f1f2ae84SBarry Smith   PetscCount         nz;
285f1f2ae84SBarry Smith   const PetscScalar *sa = NULL;
286dad3da8eSBarry Smith   PetscMPIInt        size;
28768a21331SBarry Smith   PetscInt           matproperties[4] = {0, 0, 0, 0};
288f1f2ae84SBarry Smith 
289f1f2ae84SBarry Smith   PetscFunctionBegin;
290f1f2ae84SBarry Smith   if (pc) {
291f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
292f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
293f1f2ae84SBarry Smith   }
294f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
2953ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
296f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
297dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
298dad3da8eSBarry Smith   PCMPIMatCounts[size - 1]++;
299f1f2ae84SBarry Smith   PetscCall(KSPGetOperators(ksp, NULL, &A));
300f1f2ae84SBarry Smith   PetscCall(MatMPIAIJGetNumberNonzeros(A, &nz));
301f1f2ae84SBarry Smith   PetscCall(PetscMalloc1(nz, &a));
302f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, pc ? km->NZ : NULL, pc ? km->NZdispl : NULL, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
30368a21331SBarry Smith   if (pc) {
30468a21331SBarry Smith     PetscBool isset, issymmetric, ishermitian, isspd, isstructurallysymmetric;
30568a21331SBarry Smith 
30668a21331SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
30768a21331SBarry Smith 
30868a21331SBarry Smith     PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric));
30968a21331SBarry Smith     matproperties[0] = !isset ? 0 : (issymmetric ? 1 : 2);
31068a21331SBarry Smith     PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian));
31168a21331SBarry Smith     matproperties[1] = !isset ? 0 : (ishermitian ? 1 : 2);
31268a21331SBarry Smith     PetscCall(MatIsSPDKnown(sA, &isset, &isspd));
31368a21331SBarry Smith     matproperties[2] = !isset ? 0 : (isspd ? 1 : 2);
31468a21331SBarry Smith     PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric));
31568a21331SBarry Smith     matproperties[3] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2);
31668a21331SBarry Smith   }
317f1f2ae84SBarry Smith   PetscCall(MatUpdateMPIAIJWithArray(A, a));
318f1f2ae84SBarry Smith   PetscCall(PetscFree(a));
31968a21331SBarry Smith   PetscCallMPI(MPI_Bcast(matproperties, 4, MPIU_INT, 0, comm));
32068a21331SBarry 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 */
32168a21331SBarry Smith   if (matproperties[0]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[0] == 1 ? PETSC_TRUE : PETSC_FALSE));
32268a21331SBarry Smith   if (matproperties[1]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[1] == 1 ? PETSC_TRUE : PETSC_FALSE));
32368a21331SBarry Smith   if (matproperties[2]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[2] == 1 ? PETSC_TRUE : PETSC_FALSE));
32468a21331SBarry Smith   if (matproperties[3]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE));
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326f1f2ae84SBarry Smith }
327f1f2ae84SBarry Smith 
328d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISolve(PC pc, Vec B, Vec X)
329d71ae5a4SJacob Faibussowitsch {
330f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
331f1f2ae84SBarry Smith   KSP                ksp;
332f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
333f1f2ae84SBarry Smith   const PetscScalar *sb   = NULL, *x;
334f1f2ae84SBarry Smith   PetscScalar       *b, *sx = NULL;
3355316cbedSBarry Smith   PetscInt           its, n;
3365316cbedSBarry Smith   PetscMPIInt        size;
337f1f2ae84SBarry Smith 
338f1f2ae84SBarry Smith   PetscFunctionBegin;
339f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : &ksp, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
3403ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
341f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
342f1f2ae84SBarry Smith 
343f1f2ae84SBarry Smith   /* TODO: optimize code to not require building counts/displ every time */
344f1f2ae84SBarry Smith 
345f1f2ae84SBarry Smith   /* scatterv rhs */
346dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
3475316cbedSBarry Smith   if (pc) {
3485316cbedSBarry Smith     PetscInt N;
3495316cbedSBarry Smith 
350dad3da8eSBarry Smith     PCMPISolveCounts[size - 1]++;
35168a21331SBarry Smith     PetscCall(MatGetSize(pc->pmat, &N, NULL));
3525316cbedSBarry Smith     PCMPISizes[size - 1] += N;
353f1f2ae84SBarry Smith     PetscCall(VecGetArrayRead(B, &sb));
354f1f2ae84SBarry Smith   }
355f1f2ae84SBarry Smith   PetscCall(VecGetLocalSize(ksp->vec_rhs, &n));
356f1f2ae84SBarry Smith   PetscCall(VecGetArray(ksp->vec_rhs, &b));
357f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sb, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, b, n, MPIU_SCALAR, 0, comm));
358f1f2ae84SBarry Smith   PetscCall(VecRestoreArray(ksp->vec_rhs, &b));
359f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArrayRead(B, &sb));
360f1f2ae84SBarry Smith 
3610316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
362f1f2ae84SBarry Smith   PetscCall(KSPSolve(ksp, NULL, NULL));
3630316ec64SBarry Smith   PetscCall(PetscLogStagePop());
3645316cbedSBarry Smith   PetscCall(KSPGetIterationNumber(ksp, &its));
3655316cbedSBarry Smith   PCMPIIterations[size - 1] += its;
366f1f2ae84SBarry Smith 
367f1f2ae84SBarry Smith   /* gather solution */
368f1f2ae84SBarry Smith   PetscCall(VecGetArrayRead(ksp->vec_sol, &x));
369f1f2ae84SBarry Smith   if (pc) PetscCall(VecGetArray(X, &sx));
370f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gatherv(x, n, MPIU_SCALAR, sx, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, 0, comm));
371f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArray(X, &sx));
372f1f2ae84SBarry Smith   PetscCall(VecRestoreArrayRead(ksp->vec_sol, &x));
3733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
374f1f2ae84SBarry Smith }
375f1f2ae84SBarry Smith 
376d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIDestroy(PC pc)
377d71ae5a4SJacob Faibussowitsch {
378f1f2ae84SBarry Smith   PC_MPI  *km = pc ? (PC_MPI *)pc->data : NULL;
379f1f2ae84SBarry Smith   KSP      ksp;
380f1f2ae84SBarry Smith   MPI_Comm comm = PC_MPI_COMM_WORLD;
381f1f2ae84SBarry Smith 
382f1f2ae84SBarry Smith   PetscFunctionBegin;
383f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
3843ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
385c7d372c4SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
386f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&ksp));
387c7d372c4SBarry Smith   PetscCall(PetscLogStagePop());
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
389f1f2ae84SBarry Smith }
390f1f2ae84SBarry Smith 
3913821be0aSBarry Smith PetscBool PCMPIServerActive = PETSC_FALSE;
3923821be0aSBarry Smith 
393f1f2ae84SBarry Smith /*@C
3947a99bfcaSBarry Smith   PCMPIServerBegin - starts a server that runs on the `rank != 0` MPI processes waiting to process requests for
395f1580f4eSBarry Smith   parallel `KSP` solves and management of parallel `KSP` objects.
396f1f2ae84SBarry Smith 
3973821be0aSBarry Smith   Logically Collective on all MPI processes except rank 0
398f1f2ae84SBarry Smith 
399f1580f4eSBarry Smith   Options Database Keys:
400f1f2ae84SBarry 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
4013821be0aSBarry 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
402f1f2ae84SBarry Smith 
40320f4b53cSBarry Smith   Level: developer
40420f4b53cSBarry Smith 
405f1580f4eSBarry Smith   Note:
406f1f2ae84SBarry Smith   This is normally started automatically in `PetscInitialize()` when the option is provided
407f1f2ae84SBarry Smith 
4083821be0aSBarry Smith   See `PCMPI` for information on using the solver with a `KSP` object
4093821be0aSBarry Smith 
410f1f2ae84SBarry Smith   Developer Notes:
4113821be0aSBarry Smith   When called on MPI rank 0 this sets `PETSC_COMM_WORLD` to `PETSC_COMM_SELF` to allow a main program
412f1580f4eSBarry Smith   written with `PETSC_COMM_WORLD` to run correctly on the single rank while all the ranks
413f1580f4eSBarry Smith   (that would normally be sharing `PETSC_COMM_WORLD`) to run the solver server.
414f1f2ae84SBarry Smith 
415f1580f4eSBarry Smith   Can this be integrated into the `PetscDevice` abstraction that is currently being developed?
416f1f2ae84SBarry Smith 
4173821be0aSBarry Smith   Conceivably `PCREDISTRIBUTE` could be organized in a similar manner to simplify its usage
4183821be0aSBarry Smith 
4193821be0aSBarry Smith   This could be implemented directly at the `KSP` level instead of using the `PCMPI` wrapper object
4203821be0aSBarry Smith 
421baca6076SPierre Jolivet   The code could be extended to allow an MPI + OpenMP application to use the linear solver server concept across all shared-memory
4223821be0aSBarry Smith   nodes with a single MPI process per node for the user application but multiple MPI processes per node for the linear solver.
4233821be0aSBarry Smith 
4247a99bfcaSBarry Smith   The concept could also be extended for users's callbacks for `SNES`, `TS`, and `Tao` where the `SNESSolve()` for example, runs on
4253821be0aSBarry Smith   all MPI processes but the user callback only runs on one MPI process per node.
4263821be0aSBarry Smith 
4273821be0aSBarry 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
4283821be0aSBarry Smith   the `MPI_Comm` argument from PETSc calls.
4293821be0aSBarry Smith 
4303821be0aSBarry Smith .seealso: [](sec_pcmpi), `PCMPIServerEnd()`, `PCMPI`, `KSPCheckPCMPI()`
431f1f2ae84SBarry Smith @*/
432d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerBegin(void)
433d71ae5a4SJacob Faibussowitsch {
434f1f2ae84SBarry Smith   PetscMPIInt rank;
435f1f2ae84SBarry Smith 
436f1f2ae84SBarry Smith   PetscFunctionBegin;
4379d3446b2SPierre Jolivet   PetscCall(PetscInfo(NULL, "Starting MPI Linear Solver Server\n"));
4385e1a0e3cSBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY)) {
4395e1a0e3cSBarry Smith     PetscCall(VecInitializePackage());
4405e1a0e3cSBarry Smith     PetscCall(MatInitializePackage());
4415e1a0e3cSBarry Smith     PetscCall(DMInitializePackage());
4425e1a0e3cSBarry Smith     PetscCall(PCInitializePackage());
4435e1a0e3cSBarry Smith     PetscCall(KSPInitializePackage());
4445e1a0e3cSBarry Smith     PetscCall(SNESInitializePackage());
4455e1a0e3cSBarry Smith     PetscCall(TSInitializePackage());
4465e1a0e3cSBarry Smith     PetscCall(TaoInitializePackage());
4475e1a0e3cSBarry Smith   }
448956255efSBarry Smith   PetscCall(PetscLogStageRegister("PCMPI", &PCMPIStage));
4495e1a0e3cSBarry Smith 
450f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(PC_MPI_COMM_WORLD, &rank));
451f1f2ae84SBarry Smith   if (rank == 0) {
452f1f2ae84SBarry Smith     PETSC_COMM_WORLD  = PETSC_COMM_SELF;
4533821be0aSBarry Smith     PCMPIServerActive = PETSC_TRUE;
4543ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
455f1f2ae84SBarry Smith   }
456f1f2ae84SBarry Smith 
457f1f2ae84SBarry Smith   while (PETSC_TRUE) {
458f1f2ae84SBarry Smith     PCMPICommand request = PCMPI_CREATE;
459f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
460f1f2ae84SBarry Smith     switch (request) {
461d71ae5a4SJacob Faibussowitsch     case PCMPI_CREATE:
462d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPICreate(NULL));
463d71ae5a4SJacob Faibussowitsch       break;
464d71ae5a4SJacob Faibussowitsch     case PCMPI_SET_MAT:
465d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPISetMat(NULL));
466d71ae5a4SJacob Faibussowitsch       break;
467d71ae5a4SJacob Faibussowitsch     case PCMPI_UPDATE_MAT_VALUES:
468d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPIUpdateMatValues(NULL));
469d71ae5a4SJacob Faibussowitsch       break;
470f1f2ae84SBarry Smith     case PCMPI_VIEW:
471f1f2ae84SBarry Smith       // PetscCall(PCMPIView(NULL));
472f1f2ae84SBarry Smith       break;
473d71ae5a4SJacob Faibussowitsch     case PCMPI_SOLVE:
474d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPISolve(NULL, NULL, NULL));
475d71ae5a4SJacob Faibussowitsch       break;
476d71ae5a4SJacob Faibussowitsch     case PCMPI_DESTROY:
477d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPIDestroy(NULL));
478d71ae5a4SJacob Faibussowitsch       break;
479f1f2ae84SBarry Smith     case PCMPI_EXIT:
480f1f2ae84SBarry Smith       PetscCall(PetscFinalize());
481f1f2ae84SBarry Smith       exit(0); /* not sure if this is a good idea, but cannot return because it will run users main program */
482f1f2ae84SBarry Smith       break;
483d71ae5a4SJacob Faibussowitsch     default:
484d71ae5a4SJacob Faibussowitsch       break;
485f1f2ae84SBarry Smith     }
486f1f2ae84SBarry Smith   }
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
488f1f2ae84SBarry Smith }
489f1f2ae84SBarry Smith 
490f1f2ae84SBarry Smith /*@C
491f1f2ae84SBarry Smith   PCMPIServerEnd - ends a server that runs on the rank != 0 MPI processes waiting to process requests for
492f1580f4eSBarry Smith   parallel KSP solves and management of parallel `KSP` objects.
493f1f2ae84SBarry Smith 
49420f4b53cSBarry Smith   Logically Collective on all MPI ranks except 0
49520f4b53cSBarry Smith 
49620f4b53cSBarry Smith   Level: developer
497f1f2ae84SBarry Smith 
498f1580f4eSBarry Smith   Note:
499f1f2ae84SBarry Smith   This is normally ended automatically in `PetscFinalize()` when the option is provided
500f1f2ae84SBarry Smith 
5013821be0aSBarry Smith .seealso: [](sec_pcmpi), `PCMPIServerBegin()`, `PCMPI`, `KSPCheckPCMPI()`
502f1f2ae84SBarry Smith @*/
503d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerEnd(void)
504d71ae5a4SJacob Faibussowitsch {
505f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_EXIT;
506f1f2ae84SBarry Smith 
507f1f2ae84SBarry Smith   PetscFunctionBegin;
508f1f2ae84SBarry Smith   if (PetscGlobalRank == 0) {
509f1f2ae84SBarry Smith     PetscViewer       viewer = NULL;
510f1f2ae84SBarry Smith     PetscViewerFormat format;
511f1f2ae84SBarry Smith 
512f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
513f1f2ae84SBarry Smith     PETSC_COMM_WORLD = MPI_COMM_WORLD; /* could use PC_MPI_COMM_WORLD */
514f1f2ae84SBarry Smith     PetscOptionsBegin(PETSC_COMM_SELF, NULL, "MPI linear solver server options", NULL);
515f1f2ae84SBarry Smith     PetscCall(PetscOptionsViewer("-mpi_linear_solver_server_view", "View information about system solved with the server", "PCMPI", &viewer, &format, NULL));
516f1f2ae84SBarry Smith     PetscOptionsEnd();
517f1f2ae84SBarry Smith     if (viewer) {
518f1f2ae84SBarry Smith       PetscBool isascii;
519f1f2ae84SBarry Smith 
520f1f2ae84SBarry Smith       PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
521f1f2ae84SBarry Smith       if (isascii) {
522f1f2ae84SBarry Smith         PetscMPIInt size;
5235316cbedSBarry Smith         PetscMPIInt i;
524f1f2ae84SBarry Smith 
525f1f2ae84SBarry Smith         PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
5265316cbedSBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "MPI linear solver server statistics:\n"));
5275316cbedSBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "    Ranks        KSPSolve()s     Mats        KSPs       Avg. Size      Avg. Its\n"));
5285316cbedSBarry Smith         if (PCMPIKSPCountsSeq) {
5295316cbedSBarry Smith           PetscCall(PetscViewerASCIIPrintf(viewer, "  Sequential         %" PetscInt_FMT "                         %" PetscInt_FMT "            %" PetscInt_FMT "           %" PetscInt_FMT "\n", PCMPISolveCountsSeq, PCMPIKSPCountsSeq, PCMPISizesSeq / PCMPISolveCountsSeq, PCMPIIterationsSeq / PCMPISolveCountsSeq));
530f1f2ae84SBarry Smith         }
5315316cbedSBarry Smith         for (i = 0; i < size; i++) {
5325316cbedSBarry Smith           if (PCMPIKSPCounts[i]) {
5335316cbedSBarry 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]));
5345316cbedSBarry Smith           }
5355316cbedSBarry Smith         }
536f1f2ae84SBarry Smith       }
537cd791dc2SBarry Smith       PetscCall(PetscOptionsRestoreViewer(&viewer));
538f1f2ae84SBarry Smith     }
539f1f2ae84SBarry Smith   }
540f1f2ae84SBarry Smith   PetscCall(PCMPICommsDestroy());
5413821be0aSBarry Smith   PCMPIServerActive = PETSC_FALSE;
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
543f1f2ae84SBarry Smith }
544f1f2ae84SBarry Smith 
545f1f2ae84SBarry Smith /*
546f1f2ae84SBarry Smith     This version is used in the trivial case when the MPI parallel solver server is running on just the original MPI rank 0
547f1f2ae84SBarry Smith     because, for example, the problem is small. This version is more efficient because it does not require copying any data
548f1f2ae84SBarry Smith */
549d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Seq(PC pc)
550d71ae5a4SJacob Faibussowitsch {
551f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
552f1f2ae84SBarry Smith   Mat         sA;
553f1f2ae84SBarry Smith   const char *prefix;
554f9818f3cSJose E. Roman   char       *found = NULL, *cprefix;
555f1f2ae84SBarry Smith 
556f1f2ae84SBarry Smith   PetscFunctionBegin;
557f1f2ae84SBarry Smith   PetscCall(PCGetOperators(pc, NULL, &sA));
558f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
559f1f2ae84SBarry Smith   PetscCall(KSPCreate(PETSC_COMM_SELF, &km->ksps[0]));
5603821be0aSBarry Smith   PetscCall(KSPSetNestLevel(km->ksps[0], 1));
5613821be0aSBarry Smith   PetscCall(PetscObjectSetTabLevel((PetscObject)km->ksps[0], 1));
5623821be0aSBarry Smith 
5633821be0aSBarry Smith   /* Created KSP gets prefix of PC minus the mpi_linear_solver_server_ portion */
5643821be0aSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
5653821be0aSBarry Smith   PetscCheck(prefix, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing required prefix");
5663821be0aSBarry Smith   PetscCall(PetscStrallocpy(prefix, &cprefix));
5673821be0aSBarry Smith   PetscCall(PetscStrstr(cprefix, "mpi_linear_solver_server_", &found));
5683821be0aSBarry Smith   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing mpi_linear_solver_server_ portion of prefix");
5693821be0aSBarry Smith   *found = 0;
5703821be0aSBarry Smith   PetscCall(KSPSetOptionsPrefix(km->ksps[0], cprefix));
5713821be0aSBarry Smith   PetscCall(PetscFree(cprefix));
5723821be0aSBarry Smith 
573f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(km->ksps[0], sA, sA));
574f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(km->ksps[0]));
575f1f2ae84SBarry Smith   PetscCall(KSPSetUp(km->ksps[0]));
5763ba16761SJacob Faibussowitsch   PetscCall(PetscInfo((PetscObject)pc, "MPI parallel linear solver system is being solved directly on rank 0 due to its small size\n"));
577f1f2ae84SBarry Smith   PCMPIKSPCountsSeq++;
5783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
579f1f2ae84SBarry Smith }
580f1f2ae84SBarry Smith 
581d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Seq(PC pc, Vec b, Vec x)
582d71ae5a4SJacob Faibussowitsch {
583f1f2ae84SBarry Smith   PC_MPI  *km = (PC_MPI *)pc->data;
5845316cbedSBarry Smith   PetscInt its, n;
5855316cbedSBarry Smith   Mat      A;
586f1f2ae84SBarry Smith 
587f1f2ae84SBarry Smith   PetscFunctionBegin;
588f1f2ae84SBarry Smith   PetscCall(KSPSolve(km->ksps[0], b, x));
5895316cbedSBarry Smith   PetscCall(KSPGetIterationNumber(km->ksps[0], &its));
590f1f2ae84SBarry Smith   PCMPISolveCountsSeq++;
5915316cbedSBarry Smith   PCMPIIterationsSeq += its;
5925316cbedSBarry Smith   PetscCall(KSPGetOperators(km->ksps[0], NULL, &A));
5935316cbedSBarry Smith   PetscCall(MatGetSize(A, &n, NULL));
5945316cbedSBarry Smith   PCMPISizesSeq += n;
5953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
596f1f2ae84SBarry Smith }
597f1f2ae84SBarry Smith 
598d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Seq(PC pc, PetscViewer viewer)
599d71ae5a4SJacob Faibussowitsch {
600f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
601f1f2ae84SBarry Smith 
602f1f2ae84SBarry Smith   PetscFunctionBegin;
603f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Running MPI linear solver server directly on rank 0 due to its small size\n"));
604f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
6055316cbedSBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n"));
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607f1f2ae84SBarry Smith }
608f1f2ae84SBarry Smith 
609d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Seq(PC pc)
610d71ae5a4SJacob Faibussowitsch {
611f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
612f1f2ae84SBarry Smith 
613f1f2ae84SBarry Smith   PetscFunctionBegin;
614f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&km->ksps[0]));
615f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617f1f2ae84SBarry Smith }
618f1f2ae84SBarry Smith 
619f1f2ae84SBarry Smith /*
620f1f2ae84SBarry Smith      PCSetUp_MPI - Trigger the creation of the MPI parallel PC and copy parts of the matrix and
621dd8e379bSPierre Jolivet      right-hand side to the parallel PC
622f1f2ae84SBarry Smith */
623d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_MPI(PC pc)
624d71ae5a4SJacob Faibussowitsch {
625f1f2ae84SBarry Smith   PC_MPI      *km = (PC_MPI *)pc->data;
626f1f2ae84SBarry Smith   PetscMPIInt  rank, size;
627f1f2ae84SBarry Smith   PCMPICommand request;
628f1f2ae84SBarry Smith   PetscBool    newmatrix = PETSC_FALSE;
629f1f2ae84SBarry Smith 
630f1f2ae84SBarry Smith   PetscFunctionBegin;
631f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
632f1f2ae84SBarry 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?");
633f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
634f1f2ae84SBarry Smith 
635f1f2ae84SBarry Smith   if (!pc->setupcalled) {
636f1f2ae84SBarry Smith     if (!km->alwaysuseserver) {
637f1f2ae84SBarry Smith       PetscInt n;
638f1f2ae84SBarry Smith       Mat      sA;
639f1f2ae84SBarry Smith       /* short circuit for small systems */
640f1f2ae84SBarry Smith       PetscCall(PCGetOperators(pc, &sA, &sA));
641f1f2ae84SBarry Smith       PetscCall(MatGetSize(sA, &n, NULL));
642f1f2ae84SBarry Smith       if (n < 2 * km->mincntperrank - 1 || size == 1) {
643f1f2ae84SBarry Smith         pc->ops->setup   = NULL;
644f1f2ae84SBarry Smith         pc->ops->apply   = PCApply_Seq;
645f1f2ae84SBarry Smith         pc->ops->destroy = PCDestroy_Seq;
646f1f2ae84SBarry Smith         pc->ops->view    = PCView_Seq;
647f1f2ae84SBarry Smith         PetscCall(PCSetUp_Seq(pc));
6483ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
649f1f2ae84SBarry Smith       }
650f1f2ae84SBarry Smith     }
651f1f2ae84SBarry Smith 
652f1f2ae84SBarry Smith     request = PCMPI_CREATE;
653f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
654f1f2ae84SBarry Smith     PetscCall(PCMPICreate(pc));
655f1f2ae84SBarry Smith     newmatrix = PETSC_TRUE;
6569371c9d4SSatish Balay   }
6579371c9d4SSatish Balay   if (pc->flag == DIFFERENT_NONZERO_PATTERN) newmatrix = PETSC_TRUE;
658f1f2ae84SBarry Smith 
659f1f2ae84SBarry Smith   if (newmatrix) {
6603ba16761SJacob Faibussowitsch     PetscCall(PetscInfo((PetscObject)pc, "New matrix or matrix has changed nonzero structure\n"));
661f1f2ae84SBarry Smith     request = PCMPI_SET_MAT;
662f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
663f1f2ae84SBarry Smith     PetscCall(PCMPISetMat(pc));
664f1f2ae84SBarry Smith   } else {
665bbea24aaSStefano Zampini     PetscCall(PetscInfo((PetscObject)pc, "Matrix has only changed nonzero values\n"));
666f1f2ae84SBarry Smith     request = PCMPI_UPDATE_MAT_VALUES;
667f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
668f1f2ae84SBarry Smith     PetscCall(PCMPIUpdateMatValues(pc));
669f1f2ae84SBarry Smith   }
6703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
671f1f2ae84SBarry Smith }
672f1f2ae84SBarry Smith 
673d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MPI(PC pc, Vec b, Vec x)
674d71ae5a4SJacob Faibussowitsch {
675f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_SOLVE;
676f1f2ae84SBarry Smith 
677f1f2ae84SBarry Smith   PetscFunctionBegin;
678f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
679f1f2ae84SBarry Smith   PetscCall(PCMPISolve(pc, b, x));
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
681f1f2ae84SBarry Smith }
682f1f2ae84SBarry Smith 
68366976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_MPI(PC pc)
684d71ae5a4SJacob Faibussowitsch {
685f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_DESTROY;
686f1f2ae84SBarry Smith 
687f1f2ae84SBarry Smith   PetscFunctionBegin;
688f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
689f1f2ae84SBarry Smith   PetscCall(PCMPIDestroy(pc));
690f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
6913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
692f1f2ae84SBarry Smith }
693f1f2ae84SBarry Smith 
694f1f2ae84SBarry Smith /*
695f1f2ae84SBarry Smith      PCView_MPI - Cannot call view on the MPI parallel KSP because other ranks do not have access to the viewer
696f1f2ae84SBarry Smith */
69766976f2fSJacob Faibussowitsch static PetscErrorCode PCView_MPI(PC pc, PetscViewer viewer)
698d71ae5a4SJacob Faibussowitsch {
699f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
700f1f2ae84SBarry Smith   MPI_Comm    comm;
701f1f2ae84SBarry Smith   PetscMPIInt size;
702f1f2ae84SBarry Smith 
703f1f2ae84SBarry Smith   PetscFunctionBegin;
704f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)km->ksps[0], &comm));
705f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
706f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Size of MPI communicator used for MPI parallel KSP solve %d\n", size));
707f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
70868a21331SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n"));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
710f1f2ae84SBarry Smith }
711f1f2ae84SBarry Smith 
71266976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_MPI(PC pc, PetscOptionItems *PetscOptionsObject)
713d71ae5a4SJacob Faibussowitsch {
714f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
715f1f2ae84SBarry Smith 
716f1f2ae84SBarry Smith   PetscFunctionBegin;
717f1f2ae84SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPI linear solver server options");
7183821be0aSBarry Smith   PetscCall(PetscOptionsInt("-minimum_count_per_rank", "Desired minimum number of nonzeros per rank", "None", km->mincntperrank, &km->mincntperrank, NULL));
7193821be0aSBarry 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));
720f1f2ae84SBarry Smith   PetscOptionsHeadEnd();
7213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
722f1f2ae84SBarry Smith }
723f1f2ae84SBarry Smith 
724f1f2ae84SBarry Smith /*MC
725f1580f4eSBarry Smith      PCMPI - Calls an MPI parallel `KSP` to solve a linear system from user code running on one process
726f1f2ae84SBarry Smith 
7273821be0aSBarry Smith    Options Database Keys for the Server:
728f1f2ae84SBarry 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
7293821be0aSBarry Smith -  -mpi_linear_solver_server_view - displays information about all the linear systems solved by the MPI linear solver server
730f1f2ae84SBarry Smith 
7313821be0aSBarry Smith    Options Database Keys for a specific `KSP` object
7323821be0aSBarry 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
7333821be0aSBarry 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)
7343821be0aSBarry Smith 
7353821be0aSBarry Smith    Level: developer
73620f4b53cSBarry Smith 
737f1f2ae84SBarry Smith    Notes:
73846bbbc36SPierre 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.
739f1f2ae84SBarry Smith 
740f1f2ae84SBarry Smith    It can be particularly useful for user OpenMP code or potentially user GPU code.
741f1f2ae84SBarry Smith 
742dd8e379bSPierre Jolivet    When the program is running with a single MPI process then it directly uses the provided matrix and right-hand side
7433821be0aSBarry 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.
744f1f2ae84SBarry Smith 
7453821be0aSBarry 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
7460316ec64SBarry Smith    because they are not the actual solver objects.
7470316ec64SBarry Smith 
7480316ec64SBarry 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
74968a21331SBarry 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.
7500316ec64SBarry Smith 
7513821be0aSBarry Smith    Developer Note:
7523821be0aSBarry 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
7533821be0aSBarry Smith    a `KSP`. The outer most `KSP` object is automatically set to `KSPPREONLY` and thus is not directly visible to the user.
7543821be0aSBarry Smith 
7553821be0aSBarry Smith .seealso: [](sec_pcmpi), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `PC`, `PCMPIServerBegin()`, `PCMPIServerEnd()`, `KSPCheckPCMPI()`
756f1f2ae84SBarry Smith M*/
757d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MPI(PC pc)
758d71ae5a4SJacob Faibussowitsch {
759f1f2ae84SBarry Smith   PC_MPI *km;
760f9818f3cSJose E. Roman   char   *found = NULL;
761f1f2ae84SBarry Smith 
762f1f2ae84SBarry Smith   PetscFunctionBegin;
7633821be0aSBarry Smith   PetscCall(PetscStrstr(((PetscObject)pc)->prefix, "mpi_linear_solver_server_", &found));
7643821be0aSBarry Smith   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI object prefix does not have mpi_linear_solver_server_");
7653821be0aSBarry Smith 
7663821be0aSBarry Smith   /* material from PCSetType() */
7673821be0aSBarry Smith   PetscTryTypeMethod(pc, destroy);
7683821be0aSBarry Smith   pc->ops->destroy = NULL;
7693821be0aSBarry Smith   pc->data         = NULL;
7703821be0aSBarry Smith 
7713821be0aSBarry Smith   PetscCall(PetscFunctionListDestroy(&((PetscObject)pc)->qlist));
7723821be0aSBarry Smith   PetscCall(PetscMemzero(pc->ops, sizeof(struct _PCOps)));
7733821be0aSBarry Smith   pc->modifysubmatrices  = NULL;
7743821be0aSBarry Smith   pc->modifysubmatricesP = NULL;
7753821be0aSBarry Smith   pc->setupcalled        = 0;
7763821be0aSBarry Smith 
7774dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&km));
778f1f2ae84SBarry Smith   pc->data = (void *)km;
779f1f2ae84SBarry Smith 
780f1f2ae84SBarry Smith   km->mincntperrank = 10000;
781f1f2ae84SBarry Smith 
782f1f2ae84SBarry Smith   pc->ops->setup          = PCSetUp_MPI;
783f1f2ae84SBarry Smith   pc->ops->apply          = PCApply_MPI;
784f1f2ae84SBarry Smith   pc->ops->destroy        = PCDestroy_MPI;
785f1f2ae84SBarry Smith   pc->ops->view           = PCView_MPI;
786f1f2ae84SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MPI;
7873821be0aSBarry Smith   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCMPI));
7883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
789f1f2ae84SBarry Smith }
790