xref: /petsc/src/ksp/pc/impls/mpi/pcmpi.c (revision 20f4b53cbb5e9bd9ef12b76a8697d60d197cda17)
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;
415316cbedSBarry Smith static PetscInt  PCMPIIterations[PC_MPI_MAX_RANKS], PCMPISizes[PC_MPI_MAX_RANKS], PCMPIIterationsSeq = 0, PCMPISizesSeq = 0;
42f1f2ae84SBarry Smith 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPICommsCreate(void)
44d71ae5a4SJacob Faibussowitsch {
45f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
46f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
47f1f2ae84SBarry Smith 
48f1f2ae84SBarry Smith   PetscFunctionBegin;
49f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
50f1f2ae84SBarry 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");
51f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
52f1f2ae84SBarry Smith   /* comm for size 1 is useful only for debugging */
53f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
54f1f2ae84SBarry Smith     PetscMPIInt color = rank < i + 1 ? 0 : MPI_UNDEFINED;
55f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_split(comm, color, 0, &PCMPIComms[i]));
56f1f2ae84SBarry Smith     PCMPISolveCounts[i] = 0;
57f1f2ae84SBarry Smith     PCMPIKSPCounts[i]   = 0;
585316cbedSBarry Smith     PCMPIIterations[i]  = 0;
595316cbedSBarry Smith     PCMPISizes[i]       = 0;
60f1f2ae84SBarry Smith   }
61f1f2ae84SBarry Smith   PCMPICommSet = PETSC_TRUE;
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63f1f2ae84SBarry Smith }
64f1f2ae84SBarry Smith 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPICommsDestroy(void)
66d71ae5a4SJacob Faibussowitsch {
67f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
68f1f2ae84SBarry Smith   PetscMPIInt size, rank, i;
69f1f2ae84SBarry Smith 
70f1f2ae84SBarry Smith   PetscFunctionBegin;
713ba16761SJacob Faibussowitsch   if (!PCMPICommSet) PetscFunctionReturn(PETSC_SUCCESS);
72f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
73f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
74f1f2ae84SBarry Smith   for (i = 0; i < size; i++) {
75f1f2ae84SBarry Smith     if (PCMPIComms[i] != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_free(&PCMPIComms[i]));
76f1f2ae84SBarry Smith   }
77f1f2ae84SBarry Smith   PCMPICommSet = PETSC_FALSE;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79f1f2ae84SBarry Smith }
80f1f2ae84SBarry Smith 
81d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPICreate(PC pc)
82d71ae5a4SJacob Faibussowitsch {
83f1f2ae84SBarry Smith   PC_MPI     *km   = pc ? (PC_MPI *)pc->data : NULL;
84f1f2ae84SBarry Smith   MPI_Comm    comm = PC_MPI_COMM_WORLD;
85f1f2ae84SBarry Smith   KSP         ksp;
86f1f2ae84SBarry Smith   PetscInt    N[2], mincntperrank = 0;
87f1f2ae84SBarry Smith   PetscMPIInt size;
88f1f2ae84SBarry Smith   Mat         sA;
89f1f2ae84SBarry Smith   char       *prefix;
90f1f2ae84SBarry Smith   PetscMPIInt len = 0;
91f1f2ae84SBarry Smith 
92f1f2ae84SBarry Smith   PetscFunctionBegin;
93f1f2ae84SBarry Smith   if (!PCMPICommSet) PetscCall(PCMPICommsCreate());
94f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
95f1f2ae84SBarry Smith   if (pc) {
96f1f2ae84SBarry 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"));
97f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
98f1f2ae84SBarry Smith     PetscCall(MatGetSize(sA, &N[0], &N[1]));
99f1f2ae84SBarry Smith   }
100f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(N, 2, MPIU_INT, 0, comm));
101f1f2ae84SBarry Smith 
102f1f2ae84SBarry Smith   /* choose a suitable sized MPI_Comm for the problem to be solved on */
103f1f2ae84SBarry Smith   if (km) mincntperrank = km->mincntperrank;
104f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&mincntperrank, 1, MPI_INT, 0, comm));
105f1f2ae84SBarry Smith   comm = PCMPIComms[PetscMin(size, PetscMax(1, N[0] / mincntperrank)) - 1];
106f1f2ae84SBarry Smith   if (comm == MPI_COMM_NULL) {
107f1f2ae84SBarry Smith     ksp = NULL;
1083ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
109f1f2ae84SBarry Smith   }
1100316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
111f1f2ae84SBarry Smith   PetscCall(KSPCreate(comm, &ksp));
1120316ec64SBarry Smith   PetscCall(PetscLogStagePop());
113f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gather(&ksp, 1, MPI_AINT, pc ? km->ksps : NULL, 1, MPI_AINT, 0, comm));
114f1f2ae84SBarry Smith   if (pc) {
115f1f2ae84SBarry Smith     size_t slen;
116f1f2ae84SBarry Smith 
117f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
118dad3da8eSBarry Smith     PCMPIKSPCounts[size - 1]++;
119f1f2ae84SBarry Smith     PetscCall(PCGetOptionsPrefix(pc, (const char **)&prefix));
120f1f2ae84SBarry Smith     PetscCall(PetscStrlen(prefix, &slen));
121f1f2ae84SBarry Smith     len = (PetscMPIInt)slen;
122f1f2ae84SBarry Smith   }
123f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&len, 1, MPI_INT, 0, comm));
124f1f2ae84SBarry Smith   if (len) {
12548a46eb9SPierre Jolivet     if (!pc) PetscCall(PetscMalloc1(len + 1, &prefix));
126f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(prefix, len + 1, MPI_CHAR, 0, comm));
127f1f2ae84SBarry Smith     PetscCall(KSPSetOptionsPrefix(ksp, prefix));
128f1f2ae84SBarry Smith   }
129f1f2ae84SBarry Smith   PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_"));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131f1f2ae84SBarry Smith }
132f1f2ae84SBarry Smith 
133d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISetMat(PC pc)
134d71ae5a4SJacob Faibussowitsch {
135f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
136f1f2ae84SBarry Smith   Mat                A;
13768a21331SBarry Smith   PetscInt           m, n, *ia, *ja, j, bs;
138f1f2ae84SBarry Smith   Mat                sA;
139f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
140f1f2ae84SBarry Smith   KSP                ksp;
141f1f2ae84SBarry Smith   PetscLayout        layout;
142f1f2ae84SBarry Smith   const PetscInt    *IA = NULL, *JA = NULL;
143f1f2ae84SBarry Smith   const PetscInt    *range;
144f1f2ae84SBarry Smith   PetscMPIInt       *NZ = NULL, sendcounti[PC_MPI_MAX_RANKS], displi[PC_MPI_MAX_RANKS], *NZdispl = NULL, nz, size, i;
145f1f2ae84SBarry Smith   PetscScalar       *a;
146f1f2ae84SBarry Smith   const PetscScalar *sa               = NULL;
14768a21331SBarry Smith   PetscInt           matproperties[7] = {0, 0, 0, 0, 0, 0, 0};
148f1f2ae84SBarry Smith 
149f1f2ae84SBarry Smith   PetscFunctionBegin;
150f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
1513ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
152f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
153f1f2ae84SBarry Smith   if (pc) {
15468a21331SBarry Smith     PetscBool isset, issymmetric, ishermitian, isspd, isstructurallysymmetric;
15568a21331SBarry Smith 
156f1f2ae84SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
157dad3da8eSBarry Smith     PCMPIMatCounts[size - 1]++;
158f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
15968a21331SBarry Smith     PetscCall(MatGetSize(sA, &matproperties[0], &matproperties[1]));
160dd0d27b1SBarry Smith     PetscCall(MatGetBlockSize(sA, &bs));
16168a21331SBarry Smith     matproperties[2] = bs;
16268a21331SBarry Smith     PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric));
16368a21331SBarry Smith     matproperties[3] = !isset ? 0 : (issymmetric ? 1 : 2);
16468a21331SBarry Smith     PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian));
16568a21331SBarry Smith     matproperties[4] = !isset ? 0 : (ishermitian ? 1 : 2);
16668a21331SBarry Smith     PetscCall(MatIsSPDKnown(sA, &isset, &isspd));
16768a21331SBarry Smith     matproperties[5] = !isset ? 0 : (isspd ? 1 : 2);
16868a21331SBarry Smith     PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric));
16968a21331SBarry Smith     matproperties[6] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2);
170f1f2ae84SBarry Smith   }
17168a21331SBarry Smith   PetscCallMPI(MPI_Bcast(matproperties, 7, MPIU_INT, 0, comm));
172f1f2ae84SBarry Smith 
17368a21331SBarry Smith   /* determine ownership ranges of matrix columns */
174f1f2ae84SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
17568a21331SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2]));
17668a21331SBarry Smith   PetscCall(PetscLayoutSetSize(layout, matproperties[1]));
177f1f2ae84SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
178f1f2ae84SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &n));
17968a21331SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
18068a21331SBarry Smith 
18168a21331SBarry Smith   /* determine ownership ranges of matrix rows */
18268a21331SBarry Smith   PetscCall(PetscLayoutCreate(comm, &layout));
18368a21331SBarry Smith   PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2]));
18468a21331SBarry Smith   PetscCall(PetscLayoutSetSize(layout, matproperties[0]));
18568a21331SBarry Smith   PetscCall(PetscLayoutSetUp(layout));
18668a21331SBarry Smith   PetscCall(PetscLayoutGetLocalSize(layout, &m));
187f1f2ae84SBarry Smith 
188f1f2ae84SBarry Smith   /* copy over the matrix nonzero structure and values */
189f1f2ae84SBarry Smith   if (pc) {
190f1f2ae84SBarry Smith     NZ      = km->NZ;
191f1f2ae84SBarry Smith     NZdispl = km->NZdispl;
192f1f2ae84SBarry Smith     PetscCall(PetscLayoutGetRanges(layout, &range));
193f1f2ae84SBarry Smith     PetscCall(MatGetRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
194f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
195f1f2ae84SBarry Smith       sendcounti[i] = (PetscMPIInt)(1 + range[i + 1] - range[i]);
196f1f2ae84SBarry Smith       NZ[i]         = (PetscMPIInt)(IA[range[i + 1]] - IA[range[i]]);
197f1f2ae84SBarry Smith     }
198f1f2ae84SBarry Smith     displi[0]  = 0;
199f1f2ae84SBarry Smith     NZdispl[0] = 0;
200f1f2ae84SBarry Smith     for (j = 1; j < size; j++) {
201f1f2ae84SBarry Smith       displi[j]  = displi[j - 1] + sendcounti[j - 1] - 1;
202f1f2ae84SBarry Smith       NZdispl[j] = NZdispl[j - 1] + NZ[j - 1];
203f1f2ae84SBarry Smith     }
204f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
205f1f2ae84SBarry Smith   }
206f1f2ae84SBarry Smith   PetscCall(PetscLayoutDestroy(&layout));
207f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(NZ, 1, MPI_INT, &nz, 1, MPI_INT, 0, comm));
208f1f2ae84SBarry Smith 
209f1f2ae84SBarry Smith   PetscCall(PetscMalloc3(n + 1, &ia, nz, &ja, nz, &a));
210f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(IA, sendcounti, displi, MPIU_INT, ia, n + 1, MPIU_INT, 0, comm));
211f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(JA, NZ, NZdispl, MPIU_INT, ja, nz, MPIU_INT, 0, comm));
212f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, NZ, NZdispl, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
213f1f2ae84SBarry Smith 
214f1f2ae84SBarry Smith   if (pc) {
215f1f2ae84SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
216f1f2ae84SBarry Smith     PetscCall(MatRestoreRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL));
217f1f2ae84SBarry Smith   }
218f1f2ae84SBarry Smith 
219ad540459SPierre Jolivet   for (j = 1; j < n + 1; j++) ia[j] -= ia[0];
220f1f2ae84SBarry Smith   ia[0] = 0;
2210316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
22268a21331SBarry Smith   PetscCall(MatCreateMPIAIJWithArrays(comm, m, n, matproperties[0], matproperties[1], ia, ja, a, &A));
22368a21331SBarry Smith   PetscCall(MatSetBlockSize(A, matproperties[2]));
224f1f2ae84SBarry Smith   PetscCall(MatSetOptionsPrefix(A, "mpi_"));
22568a21331SBarry Smith   if (matproperties[3]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE));
22668a21331SBarry Smith   if (matproperties[4]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[4] == 1 ? PETSC_TRUE : PETSC_FALSE));
22768a21331SBarry Smith   if (matproperties[5]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[5] == 1 ? PETSC_TRUE : PETSC_FALSE));
22868a21331SBarry Smith   if (matproperties[6]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[6] == 1 ? PETSC_TRUE : PETSC_FALSE));
229f1f2ae84SBarry Smith 
230f1f2ae84SBarry Smith   PetscCall(PetscFree3(ia, ja, a));
231f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(ksp, A, A));
232f1f2ae84SBarry Smith   if (!ksp->vec_sol) PetscCall(MatCreateVecs(A, &ksp->vec_sol, &ksp->vec_rhs));
2330316ec64SBarry Smith   PetscCall(PetscLogStagePop());
234f1f2ae84SBarry Smith   if (pc) { /* needed for scatterv/gatherv of rhs and solution */
235f1f2ae84SBarry Smith     const PetscInt *range;
236f1f2ae84SBarry Smith 
237f1f2ae84SBarry Smith     PetscCall(VecGetOwnershipRanges(ksp->vec_sol, &range));
238f1f2ae84SBarry Smith     for (i = 0; i < size; i++) {
239f1f2ae84SBarry Smith       km->sendcount[i] = (PetscMPIInt)(range[i + 1] - range[i]);
240f1f2ae84SBarry Smith       km->displ[i]     = (PetscMPIInt)range[i];
241f1f2ae84SBarry Smith     }
242f1f2ae84SBarry Smith   }
243f1f2ae84SBarry Smith   PetscCall(MatDestroy(&A));
244f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(ksp));
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
246f1f2ae84SBarry Smith }
247f1f2ae84SBarry Smith 
248d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIUpdateMatValues(PC pc)
249d71ae5a4SJacob Faibussowitsch {
250f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
251f1f2ae84SBarry Smith   KSP                ksp;
252f1f2ae84SBarry Smith   Mat                sA, A;
253f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
254f1f2ae84SBarry Smith   PetscScalar       *a;
255f1f2ae84SBarry Smith   PetscCount         nz;
256f1f2ae84SBarry Smith   const PetscScalar *sa = NULL;
257dad3da8eSBarry Smith   PetscMPIInt        size;
25868a21331SBarry Smith   PetscInt           matproperties[4] = {0, 0, 0, 0};
259f1f2ae84SBarry Smith 
260f1f2ae84SBarry Smith   PetscFunctionBegin;
261f1f2ae84SBarry Smith   if (pc) {
262f1f2ae84SBarry Smith     PetscCall(PCGetOperators(pc, &sA, &sA));
263f1f2ae84SBarry Smith     PetscCall(MatSeqAIJGetArrayRead(sA, &sa));
264f1f2ae84SBarry Smith   }
265f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
2663ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
267f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
268dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
269dad3da8eSBarry Smith   PCMPIMatCounts[size - 1]++;
270f1f2ae84SBarry Smith   PetscCall(KSPGetOperators(ksp, NULL, &A));
271f1f2ae84SBarry Smith   PetscCall(MatMPIAIJGetNumberNonzeros(A, &nz));
272f1f2ae84SBarry Smith   PetscCall(PetscMalloc1(nz, &a));
273f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sa, pc ? km->NZ : NULL, pc ? km->NZdispl : NULL, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm));
27468a21331SBarry Smith   if (pc) {
27568a21331SBarry Smith     PetscBool isset, issymmetric, ishermitian, isspd, isstructurallysymmetric;
27668a21331SBarry Smith 
27768a21331SBarry Smith     PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa));
27868a21331SBarry Smith 
27968a21331SBarry Smith     PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric));
28068a21331SBarry Smith     matproperties[0] = !isset ? 0 : (issymmetric ? 1 : 2);
28168a21331SBarry Smith     PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian));
28268a21331SBarry Smith     matproperties[1] = !isset ? 0 : (ishermitian ? 1 : 2);
28368a21331SBarry Smith     PetscCall(MatIsSPDKnown(sA, &isset, &isspd));
28468a21331SBarry Smith     matproperties[2] = !isset ? 0 : (isspd ? 1 : 2);
28568a21331SBarry Smith     PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric));
28668a21331SBarry Smith     matproperties[3] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2);
28768a21331SBarry Smith   }
288f1f2ae84SBarry Smith   PetscCall(MatUpdateMPIAIJWithArray(A, a));
289f1f2ae84SBarry Smith   PetscCall(PetscFree(a));
29068a21331SBarry Smith   PetscCallMPI(MPI_Bcast(matproperties, 4, MPIU_INT, 0, comm));
29168a21331SBarry 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 */
29268a21331SBarry Smith   if (matproperties[0]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[0] == 1 ? PETSC_TRUE : PETSC_FALSE));
29368a21331SBarry Smith   if (matproperties[1]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[1] == 1 ? PETSC_TRUE : PETSC_FALSE));
29468a21331SBarry Smith   if (matproperties[2]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[2] == 1 ? PETSC_TRUE : PETSC_FALSE));
29568a21331SBarry Smith   if (matproperties[3]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE));
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297f1f2ae84SBarry Smith }
298f1f2ae84SBarry Smith 
299d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISolve(PC pc, Vec B, Vec X)
300d71ae5a4SJacob Faibussowitsch {
301f1f2ae84SBarry Smith   PC_MPI            *km = pc ? (PC_MPI *)pc->data : NULL;
302f1f2ae84SBarry Smith   KSP                ksp;
303f1f2ae84SBarry Smith   MPI_Comm           comm = PC_MPI_COMM_WORLD;
304f1f2ae84SBarry Smith   const PetscScalar *sb   = NULL, *x;
305f1f2ae84SBarry Smith   PetscScalar       *b, *sx = NULL;
3065316cbedSBarry Smith   PetscInt           its, n;
3075316cbedSBarry Smith   PetscMPIInt        size;
308f1f2ae84SBarry Smith 
309f1f2ae84SBarry Smith   PetscFunctionBegin;
310f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : &ksp, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
3113ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
312f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
313f1f2ae84SBarry Smith 
314f1f2ae84SBarry Smith   /* TODO: optimize code to not require building counts/displ every time */
315f1f2ae84SBarry Smith 
316f1f2ae84SBarry Smith   /* scatterv rhs */
317dad3da8eSBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
3185316cbedSBarry Smith   if (pc) {
3195316cbedSBarry Smith     PetscInt N;
3205316cbedSBarry Smith 
321dad3da8eSBarry Smith     PCMPISolveCounts[size - 1]++;
32268a21331SBarry Smith     PetscCall(MatGetSize(pc->pmat, &N, NULL));
32368a21331SBarry Smith     ;
3245316cbedSBarry Smith     PCMPISizes[size - 1] += N;
325f1f2ae84SBarry Smith     PetscCall(VecGetArrayRead(B, &sb));
326f1f2ae84SBarry Smith   }
327f1f2ae84SBarry Smith   PetscCall(VecGetLocalSize(ksp->vec_rhs, &n));
328f1f2ae84SBarry Smith   PetscCall(VecGetArray(ksp->vec_rhs, &b));
329f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatterv(sb, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, b, n, MPIU_SCALAR, 0, comm));
330f1f2ae84SBarry Smith   PetscCall(VecRestoreArray(ksp->vec_rhs, &b));
331f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArrayRead(B, &sb));
332f1f2ae84SBarry Smith 
3330316ec64SBarry Smith   PetscCall(PetscLogStagePush(PCMPIStage));
334f1f2ae84SBarry Smith   PetscCall(KSPSolve(ksp, NULL, NULL));
3350316ec64SBarry Smith   PetscCall(PetscLogStagePop());
3365316cbedSBarry Smith   PetscCall(KSPGetIterationNumber(ksp, &its));
3375316cbedSBarry Smith   PCMPIIterations[size - 1] += its;
338f1f2ae84SBarry Smith 
339f1f2ae84SBarry Smith   /* gather solution */
340f1f2ae84SBarry Smith   PetscCall(VecGetArrayRead(ksp->vec_sol, &x));
341f1f2ae84SBarry Smith   if (pc) PetscCall(VecGetArray(X, &sx));
342f1f2ae84SBarry Smith   PetscCallMPI(MPI_Gatherv(x, n, MPIU_SCALAR, sx, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, 0, comm));
343f1f2ae84SBarry Smith   if (pc) PetscCall(VecRestoreArray(X, &sx));
344f1f2ae84SBarry Smith   PetscCall(VecRestoreArrayRead(ksp->vec_sol, &x));
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
346f1f2ae84SBarry Smith }
347f1f2ae84SBarry Smith 
348d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIDestroy(PC pc)
349d71ae5a4SJacob Faibussowitsch {
350f1f2ae84SBarry Smith   PC_MPI  *km = pc ? (PC_MPI *)pc->data : NULL;
351f1f2ae84SBarry Smith   KSP      ksp;
352f1f2ae84SBarry Smith   MPI_Comm comm = PC_MPI_COMM_WORLD;
353f1f2ae84SBarry Smith 
354f1f2ae84SBarry Smith   PetscFunctionBegin;
355f1f2ae84SBarry Smith   PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm));
3563ba16761SJacob Faibussowitsch   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
357f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&ksp));
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
359f1f2ae84SBarry Smith }
360f1f2ae84SBarry Smith 
361f1f2ae84SBarry Smith /*@C
362f1f2ae84SBarry Smith      PCMPIServerBegin - starts a server that runs on the rank != 0 MPI processes waiting to process requests for
363f1580f4eSBarry Smith      parallel `KSP` solves and management of parallel `KSP` objects.
364f1f2ae84SBarry Smith 
365*20f4b53cSBarry Smith      Logically Collective on all MPI ranks except 0
366f1f2ae84SBarry Smith 
367f1580f4eSBarry Smith    Options Database Keys:
368f1f2ae84SBarry 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
369f1f2ae84SBarry Smith -   -mpi_linear_solver_server_view - displays information about the linear systems solved by the MPI linear solver server
370f1f2ae84SBarry Smith 
371*20f4b53cSBarry Smith      Level: developer
372*20f4b53cSBarry Smith 
373f1580f4eSBarry Smith      Note:
374f1f2ae84SBarry Smith       This is normally started automatically in `PetscInitialize()` when the option is provided
375f1f2ae84SBarry Smith 
376f1f2ae84SBarry Smith      Developer Notes:
377f1580f4eSBarry Smith        When called on rank zero this sets `PETSC_COMM_WORLD` to `PETSC_COMM_SELF` to allow a main program
378f1580f4eSBarry Smith        written with `PETSC_COMM_WORLD` to run correctly on the single rank while all the ranks
379f1580f4eSBarry Smith        (that would normally be sharing `PETSC_COMM_WORLD`) to run the solver server.
380f1f2ae84SBarry Smith 
381f1580f4eSBarry Smith        Can this be integrated into the `PetscDevice` abstraction that is currently being developed?
382f1f2ae84SBarry Smith 
383f1580f4eSBarry Smith .seealso: `PCMPIServerEnd()`, `PCMPI`
384f1f2ae84SBarry Smith @*/
385d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerBegin(void)
386d71ae5a4SJacob Faibussowitsch {
387f1f2ae84SBarry Smith   PetscMPIInt rank;
388f1f2ae84SBarry Smith 
389f1f2ae84SBarry Smith   PetscFunctionBegin;
390f1f2ae84SBarry Smith   PetscCall(PetscInfo(NULL, "Starting MPI Linear Solver Server"));
3915e1a0e3cSBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY)) {
3925e1a0e3cSBarry Smith     PetscCall(VecInitializePackage());
3935e1a0e3cSBarry Smith     PetscCall(MatInitializePackage());
3945e1a0e3cSBarry Smith     PetscCall(DMInitializePackage());
3955e1a0e3cSBarry Smith     PetscCall(PCInitializePackage());
3965e1a0e3cSBarry Smith     PetscCall(KSPInitializePackage());
3975e1a0e3cSBarry Smith     PetscCall(SNESInitializePackage());
3985e1a0e3cSBarry Smith     PetscCall(TSInitializePackage());
3995e1a0e3cSBarry Smith     PetscCall(TaoInitializePackage());
4005e1a0e3cSBarry Smith   }
4015e1a0e3cSBarry Smith 
402f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(PC_MPI_COMM_WORLD, &rank));
403f1f2ae84SBarry Smith   if (rank == 0) {
404f1f2ae84SBarry Smith     PETSC_COMM_WORLD = PETSC_COMM_SELF;
4053ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
406f1f2ae84SBarry Smith   }
407f1f2ae84SBarry Smith 
408f1f2ae84SBarry Smith   while (PETSC_TRUE) {
409f1f2ae84SBarry Smith     PCMPICommand request = PCMPI_CREATE;
410f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
411f1f2ae84SBarry Smith     switch (request) {
412d71ae5a4SJacob Faibussowitsch     case PCMPI_CREATE:
413d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPICreate(NULL));
414d71ae5a4SJacob Faibussowitsch       break;
415d71ae5a4SJacob Faibussowitsch     case PCMPI_SET_MAT:
416d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPISetMat(NULL));
417d71ae5a4SJacob Faibussowitsch       break;
418d71ae5a4SJacob Faibussowitsch     case PCMPI_UPDATE_MAT_VALUES:
419d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPIUpdateMatValues(NULL));
420d71ae5a4SJacob Faibussowitsch       break;
421f1f2ae84SBarry Smith     case PCMPI_VIEW:
422f1f2ae84SBarry Smith       // PetscCall(PCMPIView(NULL));
423f1f2ae84SBarry Smith       break;
424d71ae5a4SJacob Faibussowitsch     case PCMPI_SOLVE:
425d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPISolve(NULL, NULL, NULL));
426d71ae5a4SJacob Faibussowitsch       break;
427d71ae5a4SJacob Faibussowitsch     case PCMPI_DESTROY:
428d71ae5a4SJacob Faibussowitsch       PetscCall(PCMPIDestroy(NULL));
429d71ae5a4SJacob Faibussowitsch       break;
430f1f2ae84SBarry Smith     case PCMPI_EXIT:
431f1f2ae84SBarry Smith       PetscCall(PetscFinalize());
432f1f2ae84SBarry Smith       exit(0); /* not sure if this is a good idea, but cannot return because it will run users main program */
433f1f2ae84SBarry Smith       break;
434d71ae5a4SJacob Faibussowitsch     default:
435d71ae5a4SJacob Faibussowitsch       break;
436f1f2ae84SBarry Smith     }
437f1f2ae84SBarry Smith   }
4383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
439f1f2ae84SBarry Smith }
440f1f2ae84SBarry Smith 
441f1f2ae84SBarry Smith /*@C
442f1f2ae84SBarry Smith      PCMPIServerEnd - ends a server that runs on the rank != 0 MPI processes waiting to process requests for
443f1580f4eSBarry Smith      parallel KSP solves and management of parallel `KSP` objects.
444f1f2ae84SBarry Smith 
445*20f4b53cSBarry Smith      Logically Collective on all MPI ranks except 0
446*20f4b53cSBarry Smith 
447*20f4b53cSBarry Smith      Level: developer
448f1f2ae84SBarry Smith 
449f1580f4eSBarry Smith      Note:
450f1f2ae84SBarry Smith       This is normally ended automatically in `PetscFinalize()` when the option is provided
451f1f2ae84SBarry Smith 
452f1580f4eSBarry Smith .seealso: `PCMPIServerBegin()`, `PCMPI`
453f1f2ae84SBarry Smith @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerEnd(void)
455d71ae5a4SJacob Faibussowitsch {
456f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_EXIT;
457f1f2ae84SBarry Smith 
458f1f2ae84SBarry Smith   PetscFunctionBegin;
459f1f2ae84SBarry Smith   if (PetscGlobalRank == 0) {
460f1f2ae84SBarry Smith     PetscViewer       viewer = NULL;
461f1f2ae84SBarry Smith     PetscViewerFormat format;
462f1f2ae84SBarry Smith 
463f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD));
464f1f2ae84SBarry Smith     PETSC_COMM_WORLD = MPI_COMM_WORLD; /* could use PC_MPI_COMM_WORLD */
465f1f2ae84SBarry Smith     PetscOptionsBegin(PETSC_COMM_SELF, NULL, "MPI linear solver server options", NULL);
466f1f2ae84SBarry Smith     PetscCall(PetscOptionsViewer("-mpi_linear_solver_server_view", "View information about system solved with the server", "PCMPI", &viewer, &format, NULL));
467f1f2ae84SBarry Smith     PetscOptionsEnd();
468f1f2ae84SBarry Smith     if (viewer) {
469f1f2ae84SBarry Smith       PetscBool isascii;
470f1f2ae84SBarry Smith 
471f1f2ae84SBarry Smith       PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
472f1f2ae84SBarry Smith       if (isascii) {
473f1f2ae84SBarry Smith         PetscMPIInt size;
4745316cbedSBarry Smith         PetscMPIInt i;
475f1f2ae84SBarry Smith 
476f1f2ae84SBarry Smith         PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
4775316cbedSBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "MPI linear solver server statistics:\n"));
4785316cbedSBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "    Ranks        KSPSolve()s     Mats        KSPs       Avg. Size      Avg. Its\n"));
4795316cbedSBarry Smith         if (PCMPIKSPCountsSeq) {
4805316cbedSBarry Smith           PetscCall(PetscViewerASCIIPrintf(viewer, "  Sequential         %" PetscInt_FMT "                         %" PetscInt_FMT "            %" PetscInt_FMT "           %" PetscInt_FMT "\n", PCMPISolveCountsSeq, PCMPIKSPCountsSeq, PCMPISizesSeq / PCMPISolveCountsSeq, PCMPIIterationsSeq / PCMPISolveCountsSeq));
481f1f2ae84SBarry Smith         }
4825316cbedSBarry Smith         for (i = 0; i < size; i++) {
4835316cbedSBarry Smith           if (PCMPIKSPCounts[i]) {
4845316cbedSBarry 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]));
4855316cbedSBarry Smith           }
4865316cbedSBarry Smith         }
487f1f2ae84SBarry Smith       }
488f1f2ae84SBarry Smith       PetscCall(PetscViewerDestroy(&viewer));
489f1f2ae84SBarry Smith     }
490f1f2ae84SBarry Smith   }
491f1f2ae84SBarry Smith   PetscCall(PCMPICommsDestroy());
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
493f1f2ae84SBarry Smith }
494f1f2ae84SBarry Smith 
495f1f2ae84SBarry Smith /*
496f1f2ae84SBarry Smith     This version is used in the trivial case when the MPI parallel solver server is running on just the original MPI rank 0
497f1f2ae84SBarry Smith     because, for example, the problem is small. This version is more efficient because it does not require copying any data
498f1f2ae84SBarry Smith */
499d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Seq(PC pc)
500d71ae5a4SJacob Faibussowitsch {
501f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
502f1f2ae84SBarry Smith   Mat         sA;
503f1f2ae84SBarry Smith   const char *prefix;
504f1f2ae84SBarry Smith 
505f1f2ae84SBarry Smith   PetscFunctionBegin;
506f1f2ae84SBarry Smith   PetscCall(PCGetOperators(pc, NULL, &sA));
507f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
508f1f2ae84SBarry Smith   PetscCall(KSPCreate(PETSC_COMM_SELF, &km->ksps[0]));
509f1f2ae84SBarry Smith   PetscCall(KSPSetOptionsPrefix(km->ksps[0], prefix));
510f1f2ae84SBarry Smith   PetscCall(KSPAppendOptionsPrefix(km->ksps[0], "mpi_"));
511f1f2ae84SBarry Smith   PetscCall(KSPSetOperators(km->ksps[0], sA, sA));
512f1f2ae84SBarry Smith   PetscCall(KSPSetFromOptions(km->ksps[0]));
513f1f2ae84SBarry Smith   PetscCall(KSPSetUp(km->ksps[0]));
5143ba16761SJacob Faibussowitsch   PetscCall(PetscInfo((PetscObject)pc, "MPI parallel linear solver system is being solved directly on rank 0 due to its small size\n"));
515f1f2ae84SBarry Smith   PCMPIKSPCountsSeq++;
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
517f1f2ae84SBarry Smith }
518f1f2ae84SBarry Smith 
519d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Seq(PC pc, Vec b, Vec x)
520d71ae5a4SJacob Faibussowitsch {
521f1f2ae84SBarry Smith   PC_MPI  *km = (PC_MPI *)pc->data;
5225316cbedSBarry Smith   PetscInt its, n;
5235316cbedSBarry Smith   Mat      A;
524f1f2ae84SBarry Smith 
525f1f2ae84SBarry Smith   PetscFunctionBegin;
526f1f2ae84SBarry Smith   PetscCall(KSPSolve(km->ksps[0], b, x));
5275316cbedSBarry Smith   PetscCall(KSPGetIterationNumber(km->ksps[0], &its));
528f1f2ae84SBarry Smith   PCMPISolveCountsSeq++;
5295316cbedSBarry Smith   PCMPIIterationsSeq += its;
5305316cbedSBarry Smith   PetscCall(KSPGetOperators(km->ksps[0], NULL, &A));
5315316cbedSBarry Smith   PetscCall(MatGetSize(A, &n, NULL));
5325316cbedSBarry Smith   PCMPISizesSeq += n;
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
534f1f2ae84SBarry Smith }
535f1f2ae84SBarry Smith 
536d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Seq(PC pc, PetscViewer viewer)
537d71ae5a4SJacob Faibussowitsch {
538f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
539f1f2ae84SBarry Smith   const char *prefix;
540f1f2ae84SBarry Smith 
541f1f2ae84SBarry Smith   PetscFunctionBegin;
542f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Running MPI linear solver server directly on rank 0 due to its small size\n"));
543f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
544f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
545f1f2ae84SBarry Smith   if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -%smpi_ksp_view to see the MPI KSP parameters ***\n", prefix));
546f1f2ae84SBarry Smith   else PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_ksp_view to see the MPI KSP parameters ***\n"));
5475316cbedSBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n"));
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549f1f2ae84SBarry Smith }
550f1f2ae84SBarry Smith 
551d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Seq(PC pc)
552d71ae5a4SJacob Faibussowitsch {
553f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
554f1f2ae84SBarry Smith 
555f1f2ae84SBarry Smith   PetscFunctionBegin;
556f1f2ae84SBarry Smith   PetscCall(KSPDestroy(&km->ksps[0]));
557f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
5583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
559f1f2ae84SBarry Smith }
560f1f2ae84SBarry Smith 
561f1f2ae84SBarry Smith /*
562f1f2ae84SBarry Smith      PCSetUp_MPI - Trigger the creation of the MPI parallel PC and copy parts of the matrix and
563f1f2ae84SBarry Smith      right hand side to the parallel PC
564f1f2ae84SBarry Smith */
565d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_MPI(PC pc)
566d71ae5a4SJacob Faibussowitsch {
567f1f2ae84SBarry Smith   PC_MPI      *km = (PC_MPI *)pc->data;
568f1f2ae84SBarry Smith   PetscMPIInt  rank, size;
569f1f2ae84SBarry Smith   PCMPICommand request;
570f1f2ae84SBarry Smith   PetscBool    newmatrix = PETSC_FALSE;
571f1f2ae84SBarry Smith 
572f1f2ae84SBarry Smith   PetscFunctionBegin;
573f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
574f1f2ae84SBarry 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?");
575f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
576f1f2ae84SBarry Smith 
577f1f2ae84SBarry Smith   if (!pc->setupcalled) {
578f1f2ae84SBarry Smith     if (!km->alwaysuseserver) {
579f1f2ae84SBarry Smith       PetscInt n;
580f1f2ae84SBarry Smith       Mat      sA;
581f1f2ae84SBarry Smith       /* short circuit for small systems */
582f1f2ae84SBarry Smith       PetscCall(PCGetOperators(pc, &sA, &sA));
583f1f2ae84SBarry Smith       PetscCall(MatGetSize(sA, &n, NULL));
584f1f2ae84SBarry Smith       if (n < 2 * km->mincntperrank - 1 || size == 1) {
585f1f2ae84SBarry Smith         pc->ops->setup   = NULL;
586f1f2ae84SBarry Smith         pc->ops->apply   = PCApply_Seq;
587f1f2ae84SBarry Smith         pc->ops->destroy = PCDestroy_Seq;
588f1f2ae84SBarry Smith         pc->ops->view    = PCView_Seq;
589f1f2ae84SBarry Smith         PetscCall(PCSetUp_Seq(pc));
5903ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
591f1f2ae84SBarry Smith       }
592f1f2ae84SBarry Smith     }
593f1f2ae84SBarry Smith 
594f1f2ae84SBarry Smith     request = PCMPI_CREATE;
595f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
596f1f2ae84SBarry Smith     PetscCall(PCMPICreate(pc));
597f1f2ae84SBarry Smith     newmatrix = PETSC_TRUE;
5989371c9d4SSatish Balay   }
5999371c9d4SSatish Balay   if (pc->flag == DIFFERENT_NONZERO_PATTERN) newmatrix = PETSC_TRUE;
600f1f2ae84SBarry Smith 
601f1f2ae84SBarry Smith   if (newmatrix) {
6023ba16761SJacob Faibussowitsch     PetscCall(PetscInfo((PetscObject)pc, "New matrix or matrix has changed nonzero structure\n"));
603f1f2ae84SBarry Smith     request = PCMPI_SET_MAT;
604f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
605f1f2ae84SBarry Smith     PetscCall(PCMPISetMat(pc));
606f1f2ae84SBarry Smith   } else {
6073ba16761SJacob Faibussowitsch     PetscCall(PetscInfo((PetscObject)pc, "Matrix has only changed nozero values\n"));
608f1f2ae84SBarry Smith     request = PCMPI_UPDATE_MAT_VALUES;
609f1f2ae84SBarry Smith     PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
610f1f2ae84SBarry Smith     PetscCall(PCMPIUpdateMatValues(pc));
611f1f2ae84SBarry Smith   }
6123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
613f1f2ae84SBarry Smith }
614f1f2ae84SBarry Smith 
615d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MPI(PC pc, Vec b, Vec x)
616d71ae5a4SJacob Faibussowitsch {
617f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_SOLVE;
618f1f2ae84SBarry Smith 
619f1f2ae84SBarry Smith   PetscFunctionBegin;
620f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
621f1f2ae84SBarry Smith   PetscCall(PCMPISolve(pc, b, x));
6223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
623f1f2ae84SBarry Smith }
624f1f2ae84SBarry Smith 
625d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_MPI(PC pc)
626d71ae5a4SJacob Faibussowitsch {
627f1f2ae84SBarry Smith   PCMPICommand request = PCMPI_DESTROY;
628f1f2ae84SBarry Smith 
629f1f2ae84SBarry Smith   PetscFunctionBegin;
630f1f2ae84SBarry Smith   PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD));
631f1f2ae84SBarry Smith   PetscCall(PCMPIDestroy(pc));
632f1f2ae84SBarry Smith   PetscCall(PetscFree(pc->data));
6333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
634f1f2ae84SBarry Smith }
635f1f2ae84SBarry Smith 
636f1f2ae84SBarry Smith /*
637f1f2ae84SBarry Smith      PCView_MPI - Cannot call view on the MPI parallel KSP because other ranks do not have access to the viewer
638f1f2ae84SBarry Smith */
639d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_MPI(PC pc, PetscViewer viewer)
640d71ae5a4SJacob Faibussowitsch {
641f1f2ae84SBarry Smith   PC_MPI     *km = (PC_MPI *)pc->data;
642f1f2ae84SBarry Smith   MPI_Comm    comm;
643f1f2ae84SBarry Smith   PetscMPIInt size;
644f1f2ae84SBarry Smith   const char *prefix;
645f1f2ae84SBarry Smith 
646f1f2ae84SBarry Smith   PetscFunctionBegin;
647f1f2ae84SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)km->ksps[0], &comm));
648f1f2ae84SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
649f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Size of MPI communicator used for MPI parallel KSP solve %d\n", size));
650f1f2ae84SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank));
651f1f2ae84SBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
652f1f2ae84SBarry Smith   if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -%smpi_ksp_view to see the MPI KSP parameters ***\n", prefix));
653f1f2ae84SBarry Smith   else PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_ksp_view to see the MPI KSP parameters ***\n"));
65468a21331SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n"));
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
656f1f2ae84SBarry Smith }
657f1f2ae84SBarry Smith 
658d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_MPI(PC pc, PetscOptionItems *PetscOptionsObject)
659d71ae5a4SJacob Faibussowitsch {
660f1f2ae84SBarry Smith   PC_MPI *km = (PC_MPI *)pc->data;
661f1f2ae84SBarry Smith 
662f1f2ae84SBarry Smith   PetscFunctionBegin;
663f1f2ae84SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPI linear solver server options");
664f1f2ae84SBarry Smith   PetscCall(PetscOptionsInt("-pc_mpi_minimum_count_per_rank", "Desired minimum number of nonzeros per rank", "None", km->mincntperrank, &km->mincntperrank, NULL));
665f1f2ae84SBarry 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));
666f1f2ae84SBarry Smith   PetscOptionsHeadEnd();
6673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
668f1f2ae84SBarry Smith }
669f1f2ae84SBarry Smith 
670f1f2ae84SBarry Smith /*MC
671f1580f4eSBarry Smith      PCMPI - Calls an MPI parallel `KSP` to solve a linear system from user code running on one process
672f1f2ae84SBarry Smith 
673f1580f4eSBarry Smith    Options Database Keys:
674f1f2ae84SBarry 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
675f1f2ae84SBarry Smith .  -mpi_linear_solver_server_view - displays information about the linear systems solved by the MPI linear solver server
676f1f2ae84SBarry Smith .  -pc_mpi_minimum_count_per_rank - sets the minimum size of the linear system per MPI rank that the solver will strive for
677f1f2ae84SBarry 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
678f1f2ae84SBarry Smith 
679*20f4b53cSBarry Smith    Level: beginner
680*20f4b53cSBarry Smith 
681f1f2ae84SBarry Smith    Notes:
682f1580f4eSBarry Smith    The options database prefix for the MPI parallel `KSP` and `PC` is -mpi_
683f1f2ae84SBarry Smith 
684f1f2ae84SBarry 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
685f1f2ae84SBarry Smith    potentially 4 to 8 MPI ranks depending on the linear system being solved, solver algorithm, and the hardware.
686f1f2ae84SBarry Smith 
687f1f2ae84SBarry Smith    It can be particularly useful for user OpenMP code or potentially user GPU code.
688f1f2ae84SBarry Smith 
689f1580f4eSBarry 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)
690f1580f4eSBarry 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.
691f1f2ae84SBarry Smith 
6920316ec64SBarry Smith    The solver options for `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
6930316ec64SBarry Smith    because they are not the actual solver objects.
6940316ec64SBarry Smith 
6950316ec64SBarry 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
69668a21331SBarry 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.
6970316ec64SBarry Smith 
698f1f2ae84SBarry Smith .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `PC`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
699f1f2ae84SBarry Smith M*/
700d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MPI(PC pc)
701d71ae5a4SJacob Faibussowitsch {
702f1f2ae84SBarry Smith   PC_MPI *km;
703f1f2ae84SBarry Smith 
704f1f2ae84SBarry Smith   PetscFunctionBegin;
7054dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&km));
706f1f2ae84SBarry Smith   pc->data = (void *)km;
707f1f2ae84SBarry Smith 
708f1f2ae84SBarry Smith   km->mincntperrank = 10000;
709f1f2ae84SBarry Smith 
710f1f2ae84SBarry Smith   pc->ops->setup          = PCSetUp_MPI;
711f1f2ae84SBarry Smith   pc->ops->apply          = PCApply_MPI;
712f1f2ae84SBarry Smith   pc->ops->destroy        = PCDestroy_MPI;
713f1f2ae84SBarry Smith   pc->ops->view           = PCView_MPI;
714f1f2ae84SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MPI;
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
716f1f2ae84SBarry Smith }
717