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