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