xref: /petsc/src/ksp/pc/impls/mpi/pcmpi.c (revision dd0d27b1b907406c548b52a6eac4547f0c73d824)
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 
5f1f2ae84SBarry Smith     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>
155e1a0e3cSBarry Smith #include <petsc.h>
16f1f2ae84SBarry Smith 
17f1f2ae84SBarry Smith #define PC_MPI_MAX_RANKS  256
18f1f2ae84SBarry Smith #define PC_MPI_COMM_WORLD MPI_COMM_WORLD
19f1f2ae84SBarry Smith 
20f1f2ae84SBarry Smith typedef struct {
21f1f2ae84SBarry Smith   KSP         ksps[PC_MPI_MAX_RANKS];                               /* The addresses of the MPI parallel KSP on each rank, NULL when not on a rank. */
22f1f2ae84SBarry Smith   PetscMPIInt sendcount[PC_MPI_MAX_RANKS], displ[PC_MPI_MAX_RANKS]; /* For scatter/gather of rhs/solution */
23f1f2ae84SBarry Smith   PetscMPIInt NZ[PC_MPI_MAX_RANKS], NZdispl[PC_MPI_MAX_RANKS];      /* For scatter of nonzero values in matrix (and nonzero column indices initially */
24f1f2ae84SBarry Smith   PetscInt    mincntperrank;                                        /* minimum number of desired nonzeros per active rank in MPI parallel KSP solve */
25f1f2ae84SBarry Smith   PetscBool   alwaysuseserver;                                      /* for debugging use the server infrastructure even if only one MPI rank is used for the solve */
26f1f2ae84SBarry Smith } PC_MPI;
27f1f2ae84SBarry Smith 
289371c9d4SSatish Balay typedef enum {
299371c9d4SSatish Balay   PCMPI_EXIT, /* exit the PC server loop, means the controlling sequential program is done */
30f1f2ae84SBarry Smith   PCMPI_CREATE,
31f1f2ae84SBarry Smith   PCMPI_SET_MAT,           /* set original matrix (or one with different nonzero pattern) */
32f1f2ae84SBarry Smith   PCMPI_UPDATE_MAT_VALUES, /* update current matrix with new nonzero values */
33f1f2ae84SBarry Smith   PCMPI_SOLVE,
34f1f2ae84SBarry Smith   PCMPI_VIEW,
35f1f2ae84SBarry Smith   PCMPI_DESTROY /* destroy a KSP that is no longer needed */
36f1f2ae84SBarry Smith } PCMPICommand;
37f1f2ae84SBarry Smith 
38f1f2ae84SBarry Smith static MPI_Comm  PCMPIComms[PC_MPI_MAX_RANKS];
39f1f2ae84SBarry Smith static PetscBool PCMPICommSet = PETSC_FALSE;
40f1f2ae84SBarry Smith static PetscInt  PCMPISolveCounts[PC_MPI_MAX_RANKS], PCMPIKSPCounts[PC_MPI_MAX_RANKS], PCMPIMatCounts[PC_MPI_MAX_RANKS], PCMPISolveCountsSeq = 0, PCMPIKSPCountsSeq = 0;
41f1f2ae84SBarry Smith 
429371c9d4SSatish Balay static PetscErrorCode PCMPICommsCreate(void) {
43f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
44f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
45f1f2ae84SBarry Smith 
46f1f2ae84SBarry Smith   PetscFunctionBegin;
47f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
48f1f2ae84SBarry 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");
49f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
50f1f2ae84SBarry Smith   /* comm for size 1 is useful only for debugging */
51f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
52f1f2ae84SBarry Smith     PetscMPIInt color = rank < i + 1 ? 0 : MPI_UNDEFINED;
53f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_split(comm, color, 0, &PCMPIComms[i]));
54f1f2ae84SBarry Smith     PCMPISolveCounts[i] = 0;
55f1f2ae84SBarry Smith     PCMPIKSPCounts[i]   = 0;
56f1f2ae84SBarry Smith   }
57f1f2ae84SBarry Smith   PCMPICommSet = PETSC_TRUE;
58f1f2ae84SBarry Smith   PetscFunctionReturn(0);
59f1f2ae84SBarry Smith }
60f1f2ae84SBarry Smith 
619371c9d4SSatish Balay PetscErrorCode PCMPICommsDestroy(void) {
62f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
63f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
64f1f2ae84SBarry Smith 
65f1f2ae84SBarry Smith   PetscFunctionBegin;
66f1f2ae84SBarry Smith   if (!PCMPICommSet) PetscFunctionReturn(0);
67f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
68f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
69f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
70f1f2ae84SBarry Smith     if (PCMPIComms[i] != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_free(&PCMPIComms[i]));
71f1f2ae84SBarry Smith   }
72f1f2ae84SBarry Smith   PCMPICommSet = PETSC_FALSE;
73f1f2ae84SBarry Smith   PetscFunctionReturn(0);
74f1f2ae84SBarry Smith }
75f1f2ae84SBarry Smith 
769371c9d4SSatish Balay static PetscErrorCode PCMPICreate(PC pc) {
77f1f2ae84SBarry Smith   PC_MPI     *km   = pc ? (PC_MPI *)pc->data : NULL;
78f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
79f1f2ae84SBarry Smith   KSP         ksp;
80f1f2ae84SBarry Smith   PetscInt    N[2], mincntperrank = 0;
81f1f2ae84SBarry Smith   PetscMPIInt size;
82f1f2ae84SBarry Smith   Mat         sA;
83f1f2ae84SBarry Smith   char       *prefix;
84f1f2ae84SBarry Smith   PetscMPIInt len = 0;
85f1f2ae84SBarry Smith 
86f1f2ae84SBarry Smith   PetscFunctionBegin;
87f1f2ae84SBarry Smith   if (!PCMPICommSet) PetscCall(PCMPICommsCreate());
88f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
89f1f2ae84SBarry Smith   if (pc) {
90f1f2ae84SBarry 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"));
91f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
92f1f2ae84SBarry Smith     PetscCall(MatGetSize(sA, &N[0], &N[1]));
93f1f2ae84SBarry Smith   }
94f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(N, 2, MPIU_INT, 0, comm));
95f1f2ae84SBarry Smith 
96f1f2ae84SBarry Smith   /* choose a suitable sized MPI_Comm for the problem to be solved on */
97f1f2ae84SBarry Smith   if (km) mincntperrank = km->mincntperrank;
98f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&mincntperrank, 1, MPI_INT, 0, comm));
99f1f2ae84SBarry Smith   comm = PCMPIComms[PetscMin(size, PetscMax(1, N[0] / mincntperrank)) - 1];
100f1f2ae84SBarry Smith   if (comm == MPI_COMM_NULL) {
101f1f2ae84SBarry Smith     ksp = NULL;
102f1f2ae84SBarry Smith     PetscFunctionReturn(0);
103f1f2ae84SBarry Smith   }
104f1f2ae84SBarry Smith   PetscCall(KSPCreate(comm, &ksp));
105f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gather(&ksp, 1, MPI_AINT, pc ? km->ksps : NULL, 1, MPI_AINT, 0, comm));
106f1f2ae84SBarry Smith   if (pc) {
107f1f2ae84SBarry Smith     size_t slen;
108f1f2ae84SBarry Smith 
109f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
110dad3da8eSBarry Smith     PCMPIKSPCounts[size - 1]++;
111f1f2ae84SBarry Smith     PetscCall(PCGetOptionsPrefix(pc, (const char **)&prefix));
112f1f2ae84SBarry Smith     PetscCall(PetscStrlen(prefix, &slen));
113f1f2ae84SBarry Smith     len = (PetscMPIInt)slen;
114f1f2ae84SBarry Smith   }
115f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&len, 1, MPI_INT, 0, comm));
116f1f2ae84SBarry Smith   if (len) {
11748a46eb9SPierre Jolivet     if (!pc) PetscCall(PetscMalloc1(len + 1, &prefix));
118f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(prefix, len + 1, MPI_CHAR, 0, comm));
119f1f2ae84SBarry Smith     PetscCall(KSPSetOptionsPrefix(ksp, prefix));
120f1f2ae84SBarry Smith   }
121f1f2ae84SBarry Smith   PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_"));
122f1f2ae84SBarry Smith   PetscFunctionReturn(0);
123f1f2ae84SBarry Smith }
124f1f2ae84SBarry Smith 
1259371c9d4SSatish Balay static PetscErrorCode PCMPISetMat(PC pc) {
126f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
127f1f2ae84SBarry Smith   Mat                A;
128*dd0d27b1SBarry Smith   PetscInt           N[2], n, *ia, *ja, j, bs;
129f1f2ae84SBarry Smith   Mat                sA;
130f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
131f1f2ae84SBarry Smith   KSP                ksp;
132f1f2ae84SBarry Smith   PetscLayout        layout;
133f1f2ae84SBarry Smith   const PetscInt    *IA = NULL, *JA = NULL;
134f1f2ae84SBarry Smith   const PetscInt    *range;
135f1f2ae84SBarry Smith   PetscMPIInt       *NZ = NULL, sendcounti[PC_MPI_MAX_RANKS], displi[PC_MPI_MAX_RANKS], *NZdispl = NULL, nz, size, i;
136f1f2ae84SBarry Smith   PetscScalar       *a;
137f1f2ae84SBarry Smith   const PetscScalar *sa = NULL;
138f1f2ae84SBarry Smith 
139f1f2ae84SBarry Smith   PetscFunctionBegin;
140f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
141f1f2ae84SBarry Smith   if (!ksp) PetscFunctionReturn(0);
142f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
143f1f2ae84SBarry Smith   if (pc) {
144f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
145dad3da8eSBarry Smith     PCMPIMatCounts[size - 1]++;
146f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
147f1f2ae84SBarry Smith     PetscCall(MatGetSize(sA, &N[0], &N[1]));
148*dd0d27b1SBarry Smith     PetscCall(MatGetBlockSize(sA, &bs));
149*dd0d27b1SBarry Smith     /* need to broadcast symmetry flags etc if set */
150f1f2ae84SBarry Smith   }
151f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(N, 2, MPIU_INT, 0, comm));
152*dd0d27b1SBarry Smith   PetscCallMPI(MPI_Bcast(&bs, 1, MPIU_INT, 0, comm));
153f1f2ae84SBarry Smith 
154f1f2ae84SBarry Smith   /* determine ownership ranges of matrix */
155f1f2ae84SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
156*dd0d27b1SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, bs));
157f1f2ae84SBarry Smith   PetscCall(PetscLayoutSetSize(layout, N[0]));
158f1f2ae84SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
159f1f2ae84SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &n));
160f1f2ae84SBarry Smith 
161f1f2ae84SBarry Smith   /* copy over the matrix nonzero structure and values */
162f1f2ae84SBarry Smith   if (pc) {
163f1f2ae84SBarry Smith     NZ      = km->NZ;
164f1f2ae84SBarry Smith     NZdispl = km->NZdispl;
165f1f2ae84SBarry Smith     PetscCall(PetscLayoutGetRanges(layout, &range));
166f1f2ae84SBarry Smith     PetscCall(MatGetRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
167f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
168f1f2ae84SBarry Smith       sendcounti[i] = (PetscMPIInt)(1 + range[i + 1] - range[i]);
169f1f2ae84SBarry Smith       NZ[i]         = (PetscMPIInt)(IA[range[i + 1]] - IA[range[i]]);
170f1f2ae84SBarry Smith     }
171f1f2ae84SBarry Smith     displi[0]  = 0;
172f1f2ae84SBarry Smith     NZdispl[0] = 0;
173f1f2ae84SBarry Smith     for (j = 1; j < size; j++) {
174f1f2ae84SBarry Smith       displi[j]  = displi[j - 1] + sendcounti[j - 1] - 1;
175f1f2ae84SBarry Smith       NZdispl[j] = NZdispl[j - 1] + NZ[j - 1];
176f1f2ae84SBarry Smith     }
177f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
178f1f2ae84SBarry Smith   }
179f1f2ae84SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
180f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(NZ, 1, MPI_INT, &nz, 1, MPI_INT, 0, comm));
181f1f2ae84SBarry Smith 
182f1f2ae84SBarry Smith   PetscCall(PetscMalloc3(n + 1, &ia, nz, &ja, nz, &a));
183f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(IA, sendcounti, displi, MPIU_INT, ia, n + 1, MPIU_INT, 0, comm));
184f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(JA, NZ, NZdispl, MPIU_INT, ja, nz, MPIU_INT, 0, comm));
185f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, NZ, NZdispl, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
186f1f2ae84SBarry Smith 
187f1f2ae84SBarry Smith   if (pc) {
188f1f2ae84SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
189f1f2ae84SBarry Smith     PetscCall(MatRestoreRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
190f1f2ae84SBarry Smith   }
191f1f2ae84SBarry Smith 
192ad540459SPierre Jolivet   for (j = 1; j < n + 1; j++) ia[j] -= ia[0];
193f1f2ae84SBarry Smith   ia[0] = 0;
194f1f2ae84SBarry Smith   PetscCall(MatCreateMPIAIJWithArrays(comm, n, n, N[0], N[0], ia, ja, a, &A));
195*dd0d27b1SBarry Smith   PetscCall(MatSetBlockSize(A, bs));
196f1f2ae84SBarry Smith   PetscCall(MatSetOptionsPrefix(A, "mpi_"));
197f1f2ae84SBarry Smith 
198f1f2ae84SBarry Smith   PetscCall(PetscFree3(ia, ja, a));
199f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(ksp, A, A));
200f1f2ae84SBarry Smith   if (!ksp->vec_sol) PetscCall(MatCreateVecs(A, &ksp->vec_sol, &ksp->vec_rhs));
201f1f2ae84SBarry Smith   if (pc) { /* needed for scatterv/gatherv of rhs and solution */
202f1f2ae84SBarry Smith     const PetscInt *range;
203f1f2ae84SBarry Smith 
204f1f2ae84SBarry Smith     PetscCall(VecGetOwnershipRanges(ksp->vec_sol, &range));
205f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
206f1f2ae84SBarry Smith       km->sendcount[i] = (PetscMPIInt)(range[i + 1] - range[i]);
207f1f2ae84SBarry Smith       km->displ[i]     = (PetscMPIInt)range[i];
208f1f2ae84SBarry Smith     }
209f1f2ae84SBarry Smith   }
210f1f2ae84SBarry Smith   PetscCall(MatDestroy(&A));
211f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(ksp));
212f1f2ae84SBarry Smith   PetscFunctionReturn(0);
213f1f2ae84SBarry Smith }
214f1f2ae84SBarry Smith 
2159371c9d4SSatish Balay static PetscErrorCode PCMPIUpdateMatValues(PC pc) {
216f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
217f1f2ae84SBarry Smith   KSP                ksp;
218f1f2ae84SBarry Smith   Mat                sA, A;
219f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
220f1f2ae84SBarry Smith   PetscScalar       *a;
221f1f2ae84SBarry Smith   PetscCount         nz;
222f1f2ae84SBarry Smith   const PetscScalar *sa = NULL;
223dad3da8eSBarry Smith   PetscMPIInt        size;
224f1f2ae84SBarry Smith 
225f1f2ae84SBarry Smith   PetscFunctionBegin;
226f1f2ae84SBarry Smith   if (pc) {
227f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
228f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
229f1f2ae84SBarry Smith   }
230f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
231f1f2ae84SBarry Smith   if (!ksp) PetscFunctionReturn(0);
232f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
233dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
234dad3da8eSBarry Smith   PCMPIMatCounts[size - 1]++;
235f1f2ae84SBarry Smith   PetscCall(KSPGetOperators(ksp, NULL, &A));
236f1f2ae84SBarry Smith   PetscCall(MatMPIAIJGetNumberNonzeros(A, &nz));
237f1f2ae84SBarry Smith   PetscCall(PetscMalloc1(nz, &a));
238f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, pc ? km->NZ : NULL, pc ? km->NZdispl : NULL, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
239f1f2ae84SBarry Smith   if (pc) PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
240f1f2ae84SBarry Smith   PetscCall(MatUpdateMPIAIJWithArray(A, a));
241f1f2ae84SBarry Smith   PetscCall(PetscFree(a));
242f1f2ae84SBarry Smith   PetscFunctionReturn(0);
243f1f2ae84SBarry Smith }
244f1f2ae84SBarry Smith 
2459371c9d4SSatish Balay static PetscErrorCode PCMPISolve(PC pc, Vec B, Vec X) {
246f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
247f1f2ae84SBarry Smith   KSP                ksp;
248f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
249f1f2ae84SBarry Smith   const PetscScalar *sb   = NULL, *x;
250f1f2ae84SBarry Smith   PetscScalar       *b, *sx = NULL;
251f1f2ae84SBarry Smith   PetscInt           n;
252f1f2ae84SBarry Smith 
253f1f2ae84SBarry Smith   PetscFunctionBegin;
254f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : &ksp, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
255f1f2ae84SBarry Smith   if (!ksp) PetscFunctionReturn(0);
256f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
257f1f2ae84SBarry Smith 
258f1f2ae84SBarry Smith   /* TODO: optimize code to not require building counts/displ everytime */
259f1f2ae84SBarry Smith 
260f1f2ae84SBarry Smith   /* scatterv rhs */
261f1f2ae84SBarry Smith   if (pc) {
262f1f2ae84SBarry Smith     PetscMPIInt size;
263f1f2ae84SBarry Smith 
264dad3da8eSBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
265dad3da8eSBarry Smith     PCMPISolveCounts[size - 1]++;
266f1f2ae84SBarry Smith     PetscCall(VecGetArrayRead(B, &sb));
267f1f2ae84SBarry Smith   }
268f1f2ae84SBarry Smith   PetscCall(VecGetLocalSize(ksp->vec_rhs, &n));
269f1f2ae84SBarry Smith   PetscCall(VecGetArray(ksp->vec_rhs, &b));
270f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sb, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, b, n, MPIU_SCALAR, 0, comm));
271f1f2ae84SBarry Smith   PetscCall(VecRestoreArray(ksp->vec_rhs, &b));
272f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArrayRead(B, &sb));
273f1f2ae84SBarry Smith 
274f1f2ae84SBarry Smith   PetscCall(KSPSolve(ksp, NULL, NULL));
275f1f2ae84SBarry Smith 
276f1f2ae84SBarry Smith   /* gather solution */
277f1f2ae84SBarry Smith   PetscCall(VecGetArrayRead(ksp->vec_sol, &x));
278f1f2ae84SBarry Smith   if (pc) PetscCall(VecGetArray(X, &sx));
279f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gatherv(x, n, MPIU_SCALAR, sx, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, 0, comm));
280f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArray(X, &sx));
281f1f2ae84SBarry Smith   PetscCall(VecRestoreArrayRead(ksp->vec_sol, &x));
282f1f2ae84SBarry Smith   PetscFunctionReturn(0);
283f1f2ae84SBarry Smith }
284f1f2ae84SBarry Smith 
2859371c9d4SSatish Balay static PetscErrorCode PCMPIDestroy(PC pc) {
286f1f2ae84SBarry Smith   PC_MPI  *km = pc ? (PC_MPI *)pc->data : NULL;
287f1f2ae84SBarry Smith   KSP      ksp;
288f1f2ae84SBarry Smith   MPI_Comm comm = PC_MPI_COMM_WORLD;
289f1f2ae84SBarry Smith 
290f1f2ae84SBarry Smith   PetscFunctionBegin;
291f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
292f1f2ae84SBarry Smith   if (!ksp) PetscFunctionReturn(0);
293f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&ksp));
294f1f2ae84SBarry Smith   PetscFunctionReturn(0);
295f1f2ae84SBarry Smith }
296f1f2ae84SBarry Smith 
297f1f2ae84SBarry Smith /*@C
298f1f2ae84SBarry Smith      PCMPIServerBegin - starts a server that runs on the rank != 0 MPI processes waiting to process requests for
299f1f2ae84SBarry Smith      parallel KSP solves and management of parallel KSP objects.
300f1f2ae84SBarry Smith 
301f1f2ae84SBarry Smith      Logically collective on all MPI ranks except 0
302f1f2ae84SBarry Smith 
303f1f2ae84SBarry Smith    Options Database:
304f1f2ae84SBarry 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
305f1f2ae84SBarry Smith -   -mpi_linear_solver_server_view - displays information about the linear systems solved by the MPI linear solver server
306f1f2ae84SBarry Smith 
307f1f2ae84SBarry Smith      Notes:
308f1f2ae84SBarry Smith       This is normally started automatically in `PetscInitialize()` when the option is provided
309f1f2ae84SBarry Smith 
310f1f2ae84SBarry Smith      Developer Notes:
311f1f2ae84SBarry Smith        When called on rank zero this sets PETSC_COMM_WORLD to PETSC_COMM_SELF to allow a main program
312f1f2ae84SBarry Smith        written with PETSC_COMM_WORLD to run correctly on the single rank while all the ranks
313f1f2ae84SBarry Smith        (that would normally be sharing PETSC_COMM_WORLD) to run the solver server.
314f1f2ae84SBarry Smith 
315f1f2ae84SBarry Smith        Can this be integrated into the PetscDevice abstraction that is currently being developed?
316f1f2ae84SBarry Smith 
317f1f2ae84SBarry Smith      Level: developer
318f1f2ae84SBarry Smith 
319f1f2ae84SBarry Smith .seealso: PCMPIServerEnd(), PCMPI
320f1f2ae84SBarry Smith @*/
3219371c9d4SSatish Balay PetscErrorCode PCMPIServerBegin(void) {
322f1f2ae84SBarry Smith   PetscMPIInt rank;
323f1f2ae84SBarry Smith 
324f1f2ae84SBarry Smith   PetscFunctionBegin;
325f1f2ae84SBarry Smith   PetscCall(PetscInfo(NULL, "Starting MPI Linear Solver Server"));
3265e1a0e3cSBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY)) {
3275e1a0e3cSBarry Smith     PetscCall(VecInitializePackage());
3285e1a0e3cSBarry Smith     PetscCall(MatInitializePackage());
3295e1a0e3cSBarry Smith     PetscCall(DMInitializePackage());
3305e1a0e3cSBarry Smith     PetscCall(PCInitializePackage());
3315e1a0e3cSBarry Smith     PetscCall(KSPInitializePackage());
3325e1a0e3cSBarry Smith     PetscCall(SNESInitializePackage());
3335e1a0e3cSBarry Smith     PetscCall(TSInitializePackage());
3345e1a0e3cSBarry Smith     PetscCall(TaoInitializePackage());
3355e1a0e3cSBarry Smith   }
3365e1a0e3cSBarry Smith 
337f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(PC_MPI_COMM_WORLD, &rank));
338f1f2ae84SBarry Smith   if (rank == 0) {
339f1f2ae84SBarry Smith     PETSC_COMM_WORLD = PETSC_COMM_SELF;
340f1f2ae84SBarry Smith     PetscFunctionReturn(0);
341f1f2ae84SBarry Smith   }
342f1f2ae84SBarry Smith 
343f1f2ae84SBarry Smith   while (PETSC_TRUE) {
344f1f2ae84SBarry Smith     PCMPICommand request = PCMPI_CREATE;
345f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
346f1f2ae84SBarry Smith     switch (request) {
3479371c9d4SSatish Balay     case PCMPI_CREATE: PetscCall(PCMPICreate(NULL)); break;
3489371c9d4SSatish Balay     case PCMPI_SET_MAT: PetscCall(PCMPISetMat(NULL)); break;
3499371c9d4SSatish Balay     case PCMPI_UPDATE_MAT_VALUES: PetscCall(PCMPIUpdateMatValues(NULL)); break;
350f1f2ae84SBarry Smith     case PCMPI_VIEW:
351f1f2ae84SBarry Smith       // PetscCall(PCMPIView(NULL));
352f1f2ae84SBarry Smith       break;
3539371c9d4SSatish Balay     case PCMPI_SOLVE: PetscCall(PCMPISolve(NULL, NULL, NULL)); break;
3549371c9d4SSatish Balay     case PCMPI_DESTROY: PetscCall(PCMPIDestroy(NULL)); break;
355f1f2ae84SBarry Smith     case PCMPI_EXIT:
356f1f2ae84SBarry Smith       PetscCall(PetscFinalize());
357f1f2ae84SBarry Smith       exit(0); /* not sure if this is a good idea, but cannot return because it will run users main program */
358f1f2ae84SBarry Smith       break;
3599371c9d4SSatish Balay     default: break;
360f1f2ae84SBarry Smith     }
361f1f2ae84SBarry Smith   }
362f1f2ae84SBarry Smith   PetscFunctionReturn(0);
363f1f2ae84SBarry Smith }
364f1f2ae84SBarry Smith 
365f1f2ae84SBarry Smith /*@C
366f1f2ae84SBarry Smith      PCMPIServerEnd - ends a server that runs on the rank != 0 MPI processes waiting to process requests for
367f1f2ae84SBarry Smith      parallel KSP solves and management of parallel KSP objects.
368f1f2ae84SBarry Smith 
369f1f2ae84SBarry Smith      Logically collective on all MPI ranks except 0
370f1f2ae84SBarry Smith 
371f1f2ae84SBarry Smith      Notes:
372f1f2ae84SBarry Smith       This is normally ended automatically in `PetscFinalize()` when the option is provided
373f1f2ae84SBarry Smith 
374f1f2ae84SBarry Smith      Level: developer
375f1f2ae84SBarry Smith 
376f1f2ae84SBarry Smith .seealso: PCMPIServerBegin(), PCMPI
377f1f2ae84SBarry Smith @*/
3789371c9d4SSatish Balay PetscErrorCode PCMPIServerEnd(void) {
379f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_EXIT;
380f1f2ae84SBarry Smith 
381f1f2ae84SBarry Smith   PetscFunctionBegin;
382f1f2ae84SBarry Smith   if (PetscGlobalRank == 0) {
383f1f2ae84SBarry Smith     PetscViewer       viewer = NULL;
384f1f2ae84SBarry Smith     PetscViewerFormat format;
385f1f2ae84SBarry Smith 
386f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
387f1f2ae84SBarry Smith     PETSC_COMM_WORLD = MPI_COMM_WORLD; /* could use PC_MPI_COMM_WORLD */
388f1f2ae84SBarry Smith     PetscOptionsBegin(PETSC_COMM_SELF, NULL, "MPI linear solver server options", NULL);
389f1f2ae84SBarry Smith     PetscCall(PetscOptionsViewer("-mpi_linear_solver_server_view", "View information about system solved with the server", "PCMPI", &viewer, &format, NULL));
390f1f2ae84SBarry Smith     PetscOptionsEnd();
391f1f2ae84SBarry Smith     if (viewer) {
392f1f2ae84SBarry Smith       PetscBool isascii;
393f1f2ae84SBarry Smith 
394f1f2ae84SBarry Smith       PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
395f1f2ae84SBarry Smith       if (isascii) {
396f1f2ae84SBarry Smith         PetscMPIInt size;
397f1f2ae84SBarry Smith         PetscInt    i, ksp = 0, mat = 0, solve = 0;
398f1f2ae84SBarry Smith 
399f1f2ae84SBarry Smith         PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
400dad3da8eSBarry Smith         for (i = 0; i < size; i++) {
401f1f2ae84SBarry Smith           ksp += PCMPIKSPCounts[i];
402f1f2ae84SBarry Smith           mat += PCMPIMatCounts[i];
403f1f2ae84SBarry Smith           solve += PCMPISolveCounts[i];
404f1f2ae84SBarry Smith         }
405f1f2ae84SBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "MPI linear solver server:\n"));
406f1f2ae84SBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "  Parallel KSPs  %" PetscInt_FMT " Mats  %" PetscInt_FMT " Solves %" PetscInt_FMT "\n", ksp, mat, solve));
407f1f2ae84SBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "  Sequential KSPs  %" PetscInt_FMT " Solves %" PetscInt_FMT "\n", PCMPIKSPCountsSeq, PCMPISolveCountsSeq));
408f1f2ae84SBarry Smith       }
409f1f2ae84SBarry Smith       PetscCall(PetscViewerDestroy(&viewer));
410f1f2ae84SBarry Smith     }
411f1f2ae84SBarry Smith   }
412f1f2ae84SBarry Smith   PetscCall(PCMPICommsDestroy());
413f1f2ae84SBarry Smith   PetscFunctionReturn(0);
414f1f2ae84SBarry Smith }
415f1f2ae84SBarry Smith 
416f1f2ae84SBarry Smith /*
417f1f2ae84SBarry Smith     This version is used in the trivial case when the MPI parallel solver server is running on just the original MPI rank 0
418f1f2ae84SBarry Smith     because, for example, the problem is small. This version is more efficient because it does not require copying any data
419f1f2ae84SBarry Smith */
4209371c9d4SSatish Balay static PetscErrorCode PCSetUp_Seq(PC pc) {
421f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
422f1f2ae84SBarry Smith   Mat         sA;
423f1f2ae84SBarry Smith   const char *prefix;
424f1f2ae84SBarry Smith 
425f1f2ae84SBarry Smith   PetscFunctionBegin;
426f1f2ae84SBarry Smith   PetscCall(PCGetOperators(pc, NULL, &sA));
427f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
428f1f2ae84SBarry Smith   PetscCall(KSPCreate(PETSC_COMM_SELF, &km->ksps[0]));
429f1f2ae84SBarry Smith   PetscCall(KSPSetOptionsPrefix(km->ksps[0], prefix));
430f1f2ae84SBarry Smith   PetscCall(KSPAppendOptionsPrefix(km->ksps[0], "mpi_"));
431f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(km->ksps[0], sA, sA));
432f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(km->ksps[0]));
433f1f2ae84SBarry Smith   PetscCall(KSPSetUp(km->ksps[0]));
434f1f2ae84SBarry Smith   PetscInfo((PetscObject)pc, "MPI parallel linear solver system is being solved directly on rank 0 due to its small size\n");
435f1f2ae84SBarry Smith   PCMPIKSPCountsSeq++;
436f1f2ae84SBarry Smith   PetscFunctionReturn(0);
437f1f2ae84SBarry Smith }
438f1f2ae84SBarry Smith 
4399371c9d4SSatish Balay static PetscErrorCode PCApply_Seq(PC pc, Vec b, Vec x) {
440f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
441f1f2ae84SBarry Smith 
442f1f2ae84SBarry Smith   PetscFunctionBegin;
443f1f2ae84SBarry Smith   PetscCall(KSPSolve(km->ksps[0], b, x));
444f1f2ae84SBarry Smith   PCMPISolveCountsSeq++;
445f1f2ae84SBarry Smith   PetscFunctionReturn(0);
446f1f2ae84SBarry Smith }
447f1f2ae84SBarry Smith 
4489371c9d4SSatish Balay static PetscErrorCode PCView_Seq(PC pc, PetscViewer viewer) {
449f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
450f1f2ae84SBarry Smith   const char *prefix;
451f1f2ae84SBarry Smith 
452f1f2ae84SBarry Smith   PetscFunctionBegin;
453f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Running MPI linear solver server directly on rank 0 due to its small size\n"));
454f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
455f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
456f1f2ae84SBarry Smith   if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -%smpi_ksp_view to see the MPI KSP parameters***\n", prefix));
457f1f2ae84SBarry Smith   else PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_ksp_view to see the MPI KSP parameters***\n"));
458f1f2ae84SBarry Smith   PetscFunctionReturn(0);
459f1f2ae84SBarry Smith }
460f1f2ae84SBarry Smith 
4619371c9d4SSatish Balay static PetscErrorCode PCDestroy_Seq(PC pc) {
462f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
463f1f2ae84SBarry Smith 
464f1f2ae84SBarry Smith   PetscFunctionBegin;
465f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&km->ksps[0]));
466f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
467f1f2ae84SBarry Smith   PetscFunctionReturn(0);
468f1f2ae84SBarry Smith }
469f1f2ae84SBarry Smith 
470f1f2ae84SBarry Smith /*
471f1f2ae84SBarry Smith      PCSetUp_MPI - Trigger the creation of the MPI parallel PC and copy parts of the matrix and
472f1f2ae84SBarry Smith      right hand side to the parallel PC
473f1f2ae84SBarry Smith */
4749371c9d4SSatish Balay static PetscErrorCode PCSetUp_MPI(PC pc) {
475f1f2ae84SBarry Smith   PC_MPI      *km = (PC_MPI *)pc->data;
476f1f2ae84SBarry Smith   PetscMPIInt  rank, size;
477f1f2ae84SBarry Smith   PCMPICommand request;
478f1f2ae84SBarry Smith   PetscBool    newmatrix = PETSC_FALSE;
479f1f2ae84SBarry Smith 
480f1f2ae84SBarry Smith   PetscFunctionBegin;
481f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
482f1f2ae84SBarry 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?");
483f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
484f1f2ae84SBarry Smith 
485f1f2ae84SBarry Smith   if (!pc->setupcalled) {
486f1f2ae84SBarry Smith     if (!km->alwaysuseserver) {
487f1f2ae84SBarry Smith       PetscInt n;
488f1f2ae84SBarry Smith       Mat      sA;
489f1f2ae84SBarry Smith       /* short circuit for small systems */
490f1f2ae84SBarry Smith       PetscCall(PCGetOperators(pc, &sA, &sA));
491f1f2ae84SBarry Smith       PetscCall(MatGetSize(sA, &n, NULL));
492f1f2ae84SBarry Smith       if (n < 2 * km->mincntperrank - 1 || size == 1) {
493f1f2ae84SBarry Smith         pc->ops->setup   = NULL;
494f1f2ae84SBarry Smith         pc->ops->apply   = PCApply_Seq;
495f1f2ae84SBarry Smith         pc->ops->destroy = PCDestroy_Seq;
496f1f2ae84SBarry Smith         pc->ops->view    = PCView_Seq;
497f1f2ae84SBarry Smith         PetscCall(PCSetUp_Seq(pc));
498f1f2ae84SBarry Smith         PetscFunctionReturn(0);
499f1f2ae84SBarry Smith       }
500f1f2ae84SBarry Smith     }
501f1f2ae84SBarry Smith 
502f1f2ae84SBarry Smith     request = PCMPI_CREATE;
503f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
504f1f2ae84SBarry Smith     PetscCall(PCMPICreate(pc));
505f1f2ae84SBarry Smith     newmatrix = PETSC_TRUE;
5069371c9d4SSatish Balay   }
5079371c9d4SSatish Balay   if (pc->flag == DIFFERENT_NONZERO_PATTERN) newmatrix = PETSC_TRUE;
508f1f2ae84SBarry Smith 
509f1f2ae84SBarry Smith   if (newmatrix) {
510f1f2ae84SBarry Smith     PetscInfo((PetscObject)pc, "New matrix or matrix has changed nonzero structure\n");
511f1f2ae84SBarry Smith     request = PCMPI_SET_MAT;
512f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
513f1f2ae84SBarry Smith     PetscCall(PCMPISetMat(pc));
514f1f2ae84SBarry Smith   } else {
515f1f2ae84SBarry Smith     PetscInfo((PetscObject)pc, "Matrix has only changed nozero values\n");
516f1f2ae84SBarry Smith     request = PCMPI_UPDATE_MAT_VALUES;
517f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
518f1f2ae84SBarry Smith     PetscCall(PCMPIUpdateMatValues(pc));
519f1f2ae84SBarry Smith   }
520f1f2ae84SBarry Smith   PetscFunctionReturn(0);
521f1f2ae84SBarry Smith }
522f1f2ae84SBarry Smith 
5239371c9d4SSatish Balay static PetscErrorCode PCApply_MPI(PC pc, Vec b, Vec x) {
524f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_SOLVE;
525f1f2ae84SBarry Smith 
526f1f2ae84SBarry Smith   PetscFunctionBegin;
527f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
528f1f2ae84SBarry Smith   PetscCall(PCMPISolve(pc, b, x));
529f1f2ae84SBarry Smith   PetscFunctionReturn(0);
530f1f2ae84SBarry Smith }
531f1f2ae84SBarry Smith 
5329371c9d4SSatish Balay PetscErrorCode PCDestroy_MPI(PC pc) {
533f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_DESTROY;
534f1f2ae84SBarry Smith 
535f1f2ae84SBarry Smith   PetscFunctionBegin;
536f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
537f1f2ae84SBarry Smith   PetscCall(PCMPIDestroy(pc));
538f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
539f1f2ae84SBarry Smith   PetscFunctionReturn(0);
540f1f2ae84SBarry Smith }
541f1f2ae84SBarry Smith 
542f1f2ae84SBarry Smith /*
543f1f2ae84SBarry Smith      PCView_MPI - Cannot call view on the MPI parallel KSP because other ranks do not have access to the viewer
544f1f2ae84SBarry Smith */
5459371c9d4SSatish Balay PetscErrorCode PCView_MPI(PC pc, PetscViewer viewer) {
546f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
547f1f2ae84SBarry Smith   MPI_Comm    comm;
548f1f2ae84SBarry Smith   PetscMPIInt size;
549f1f2ae84SBarry Smith   const char *prefix;
550f1f2ae84SBarry Smith 
551f1f2ae84SBarry Smith   PetscFunctionBegin;
552f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)km->ksps[0], &comm));
553f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
554f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Size of MPI communicator used for MPI parallel KSP solve %d\n", size));
555f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
556f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
557f1f2ae84SBarry Smith   if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -%smpi_ksp_view to see the MPI KSP parameters***\n", prefix));
558f1f2ae84SBarry Smith   else PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_ksp_view to see the MPI KSP parameters***\n"));
559f1f2ae84SBarry Smith   PetscFunctionReturn(0);
560f1f2ae84SBarry Smith }
561f1f2ae84SBarry Smith 
5629371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_MPI(PC pc, PetscOptionItems *PetscOptionsObject) {
563f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
564f1f2ae84SBarry Smith 
565f1f2ae84SBarry Smith   PetscFunctionBegin;
566f1f2ae84SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPI linear solver server options");
567f1f2ae84SBarry Smith   PetscCall(PetscOptionsInt("-pc_mpi_minimum_count_per_rank", "Desired minimum number of nonzeros per rank", "None", km->mincntperrank, &km->mincntperrank, NULL));
568f1f2ae84SBarry Smith   PetscCall(PetscOptionsBool("-pc_mpi_always_use_server", "Use the server even if only one rank is used for the solve (for debugging)", "None", km->alwaysuseserver, &km->alwaysuseserver, NULL));
569f1f2ae84SBarry Smith   PetscOptionsHeadEnd();
570f1f2ae84SBarry Smith   PetscFunctionReturn(0);
571f1f2ae84SBarry Smith }
572f1f2ae84SBarry Smith 
573f1f2ae84SBarry Smith /*MC
574f1f2ae84SBarry Smith      PCMPI - Calls an MPI parallel KSP to solve a linear system from user code running on one process
575f1f2ae84SBarry Smith 
576f1f2ae84SBarry Smith    Level: beginner
577f1f2ae84SBarry Smith 
578f1f2ae84SBarry Smith    Options Database:
579f1f2ae84SBarry 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
580f1f2ae84SBarry Smith .  -mpi_linear_solver_server_view - displays information about the linear systems solved by the MPI linear solver server
581f1f2ae84SBarry Smith .  -pc_mpi_minimum_count_per_rank - sets the minimum size of the linear system per MPI rank that the solver will strive for
582f1f2ae84SBarry Smith -  -pc_mpi_always_use_server - use the server solver code even if the particular system is only solved on the process, this option is only for debugging and testing purposes
583f1f2ae84SBarry Smith 
584f1f2ae84SBarry Smith    Notes:
585f1f2ae84SBarry Smith    The options database prefix for the MPI parallel KSP and PC is -mpi_
586f1f2ae84SBarry Smith 
587f1f2ae84SBarry Smith    The MPI linear solver server will not support scaling user code to utilize extremely large numbers of MPI ranks but should give reasonable speedup for
588f1f2ae84SBarry Smith    potentially 4 to 8 MPI ranks depending on the linear system being solved, solver algorithm, and the hardware.
589f1f2ae84SBarry Smith 
590f1f2ae84SBarry Smith    It can be particularly useful for user OpenMP code or potentially user GPU code.
591f1f2ae84SBarry Smith 
592f1f2ae84SBarry Smith    When the program is running with a single MPI rank then this directly uses the provided matrix and right hand side (still in a KSP with the options prefix of -mpi)
593f1f2ae84SBarry Smith    and does not need to distribute the matrix and vector to the various MPI ranks; thus it incurs no extra overhead over just using the KSP directly.
594f1f2ae84SBarry Smith 
595f1f2ae84SBarry Smith .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `PC`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
596f1f2ae84SBarry Smith 
597f1f2ae84SBarry Smith M*/
5989371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_MPI(PC pc) {
599f1f2ae84SBarry Smith   PC_MPI *km;
600f1f2ae84SBarry Smith 
601f1f2ae84SBarry Smith   PetscFunctionBegin;
602f1f2ae84SBarry Smith   PetscCall(PetscNewLog(pc, &km));
603f1f2ae84SBarry Smith   pc->data = (void *)km;
604f1f2ae84SBarry Smith 
605f1f2ae84SBarry Smith   km->mincntperrank = 10000;
606f1f2ae84SBarry Smith 
607f1f2ae84SBarry Smith   pc->ops->setup          = PCSetUp_MPI;
608f1f2ae84SBarry Smith   pc->ops->apply          = PCApply_MPI;
609f1f2ae84SBarry Smith   pc->ops->destroy        = PCDestroy_MPI;
610f1f2ae84SBarry Smith   pc->ops->view           = PCView_MPI;
611f1f2ae84SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MPI;
612f1f2ae84SBarry Smith   PetscFunctionReturn(0);
613f1f2ae84SBarry Smith }
614