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