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; 893821be0aSBarry Smith char *cprefix = NULL; 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)); 1123821be0aSBarry Smith PetscCall(KSPSetNestLevel(ksp, 1)); 1133821be0aSBarry Smith PetscCall(PetscObjectSetTabLevel((PetscObject)ksp, 1)); 1140316ec64SBarry Smith PetscCall(PetscLogStagePop()); 115f1f2ae84SBarry Smith PetscCallMPI(MPI_Gather(&ksp, 1, MPI_AINT, pc ? km->ksps : NULL, 1, MPI_AINT, 0, comm)); 116f1f2ae84SBarry Smith if (pc) { 117f1f2ae84SBarry Smith size_t slen; 1183821be0aSBarry Smith const char *prefix = NULL; 119f9818f3cSJose E. Roman char *found = NULL; 120f1f2ae84SBarry Smith 121f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 122dad3da8eSBarry Smith PCMPIKSPCounts[size - 1]++; 1233821be0aSBarry Smith /* Created KSP gets prefix of PC minus the mpi_linear_solver_server_ portion */ 1243821be0aSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1253821be0aSBarry Smith PetscCheck(prefix, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing required prefix"); 1263821be0aSBarry Smith PetscCall(PetscStrallocpy(prefix, &cprefix)); 1273821be0aSBarry Smith PetscCall(PetscStrstr(cprefix, "mpi_linear_solver_server_", &found)); 1283821be0aSBarry Smith PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing mpi_linear_solver_server_ portion of prefix"); 1293821be0aSBarry Smith *found = 0; 1303821be0aSBarry Smith PetscCall(PetscStrlen(cprefix, &slen)); 131f1f2ae84SBarry Smith len = (PetscMPIInt)slen; 132f1f2ae84SBarry Smith } 133f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&len, 1, MPI_INT, 0, comm)); 134f1f2ae84SBarry Smith if (len) { 1353821be0aSBarry Smith if (!pc) PetscCall(PetscMalloc1(len + 1, &cprefix)); 1363821be0aSBarry Smith PetscCallMPI(MPI_Bcast(cprefix, len + 1, MPI_CHAR, 0, comm)); 1373821be0aSBarry Smith PetscCall(KSPSetOptionsPrefix(ksp, cprefix)); 138f1f2ae84SBarry Smith } 1393821be0aSBarry Smith PetscCall(PetscFree(cprefix)); 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141f1f2ae84SBarry Smith } 142f1f2ae84SBarry Smith 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISetMat(PC pc) 144d71ae5a4SJacob Faibussowitsch { 145f1f2ae84SBarry Smith PC_MPI *km = pc ? (PC_MPI *)pc->data : NULL; 146f1f2ae84SBarry Smith Mat A; 14768a21331SBarry Smith PetscInt m, n, *ia, *ja, j, bs; 148f1f2ae84SBarry Smith Mat sA; 149f1f2ae84SBarry Smith MPI_Comm comm = PC_MPI_COMM_WORLD; 150f1f2ae84SBarry Smith KSP ksp; 151f1f2ae84SBarry Smith PetscLayout layout; 152f1f2ae84SBarry Smith const PetscInt *IA = NULL, *JA = NULL; 153f1f2ae84SBarry Smith const PetscInt *range; 154f1f2ae84SBarry Smith PetscMPIInt *NZ = NULL, sendcounti[PC_MPI_MAX_RANKS], displi[PC_MPI_MAX_RANKS], *NZdispl = NULL, nz, size, i; 155f1f2ae84SBarry Smith PetscScalar *a; 156f1f2ae84SBarry Smith const PetscScalar *sa = NULL; 1573821be0aSBarry Smith PetscInt matproperties[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 1583821be0aSBarry Smith char *cprefix; 159f1f2ae84SBarry Smith 160f1f2ae84SBarry Smith PetscFunctionBegin; 161f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm)); 1623ba16761SJacob Faibussowitsch if (!ksp) PetscFunctionReturn(PETSC_SUCCESS); 163f1f2ae84SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm)); 164f1f2ae84SBarry Smith if (pc) { 16568a21331SBarry Smith PetscBool isset, issymmetric, ishermitian, isspd, isstructurallysymmetric; 1663821be0aSBarry Smith const char *prefix; 1673821be0aSBarry Smith size_t clen; 16868a21331SBarry Smith 169f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 170dad3da8eSBarry Smith PCMPIMatCounts[size - 1]++; 171f1f2ae84SBarry Smith PetscCall(PCGetOperators(pc, &sA, &sA)); 17268a21331SBarry Smith PetscCall(MatGetSize(sA, &matproperties[0], &matproperties[1])); 173dd0d27b1SBarry Smith PetscCall(MatGetBlockSize(sA, &bs)); 17468a21331SBarry Smith matproperties[2] = bs; 17568a21331SBarry Smith PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric)); 17668a21331SBarry Smith matproperties[3] = !isset ? 0 : (issymmetric ? 1 : 2); 17768a21331SBarry Smith PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian)); 17868a21331SBarry Smith matproperties[4] = !isset ? 0 : (ishermitian ? 1 : 2); 17968a21331SBarry Smith PetscCall(MatIsSPDKnown(sA, &isset, &isspd)); 18068a21331SBarry Smith matproperties[5] = !isset ? 0 : (isspd ? 1 : 2); 18168a21331SBarry Smith PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric)); 18268a21331SBarry Smith matproperties[6] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2); 1837a99bfcaSBarry Smith /* Created Mat gets prefix of input Mat PLUS the mpi_linear_solver_server_ portion */ 1843821be0aSBarry Smith PetscCall(MatGetOptionsPrefix(sA, &prefix)); 1853821be0aSBarry Smith PetscCall(PetscStrallocpy(prefix, &cprefix)); 1863821be0aSBarry Smith PetscCall(PetscStrlen(cprefix, &clen)); 1873821be0aSBarry Smith matproperties[7] = (PetscInt)clen; 188f1f2ae84SBarry Smith } 1893821be0aSBarry Smith PetscCallMPI(MPI_Bcast(matproperties, PETSC_STATIC_ARRAY_LENGTH(matproperties), MPIU_INT, 0, comm)); 190f1f2ae84SBarry Smith 19168a21331SBarry Smith /* determine ownership ranges of matrix columns */ 192f1f2ae84SBarry Smith PetscCall(PetscLayoutCreate(comm, &layout)); 19368a21331SBarry Smith PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2])); 19468a21331SBarry Smith PetscCall(PetscLayoutSetSize(layout, matproperties[1])); 195f1f2ae84SBarry Smith PetscCall(PetscLayoutSetUp(layout)); 196f1f2ae84SBarry Smith PetscCall(PetscLayoutGetLocalSize(layout, &n)); 19768a21331SBarry Smith PetscCall(PetscLayoutDestroy(&layout)); 19868a21331SBarry Smith 19968a21331SBarry Smith /* determine ownership ranges of matrix rows */ 20068a21331SBarry Smith PetscCall(PetscLayoutCreate(comm, &layout)); 20168a21331SBarry Smith PetscCall(PetscLayoutSetBlockSize(layout, matproperties[2])); 20268a21331SBarry Smith PetscCall(PetscLayoutSetSize(layout, matproperties[0])); 20368a21331SBarry Smith PetscCall(PetscLayoutSetUp(layout)); 20468a21331SBarry Smith PetscCall(PetscLayoutGetLocalSize(layout, &m)); 205f1f2ae84SBarry Smith 206f1f2ae84SBarry Smith /* copy over the matrix nonzero structure and values */ 207f1f2ae84SBarry Smith if (pc) { 208f1f2ae84SBarry Smith NZ = km->NZ; 209f1f2ae84SBarry Smith NZdispl = km->NZdispl; 210f1f2ae84SBarry Smith PetscCall(PetscLayoutGetRanges(layout, &range)); 211f1f2ae84SBarry Smith PetscCall(MatGetRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL)); 212f1f2ae84SBarry Smith for (i = 0; i < size; i++) { 213f1f2ae84SBarry Smith sendcounti[i] = (PetscMPIInt)(1 + range[i + 1] - range[i]); 214f1f2ae84SBarry Smith NZ[i] = (PetscMPIInt)(IA[range[i + 1]] - IA[range[i]]); 215f1f2ae84SBarry Smith } 216f1f2ae84SBarry Smith displi[0] = 0; 217f1f2ae84SBarry Smith NZdispl[0] = 0; 218f1f2ae84SBarry Smith for (j = 1; j < size; j++) { 219f1f2ae84SBarry Smith displi[j] = displi[j - 1] + sendcounti[j - 1] - 1; 220f1f2ae84SBarry Smith NZdispl[j] = NZdispl[j - 1] + NZ[j - 1]; 221f1f2ae84SBarry Smith } 222f1f2ae84SBarry Smith PetscCall(MatSeqAIJGetArrayRead(sA, &sa)); 223f1f2ae84SBarry Smith } 224f1f2ae84SBarry Smith PetscCall(PetscLayoutDestroy(&layout)); 225f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatter(NZ, 1, MPI_INT, &nz, 1, MPI_INT, 0, comm)); 226f1f2ae84SBarry Smith 227f1f2ae84SBarry Smith PetscCall(PetscMalloc3(n + 1, &ia, nz, &ja, nz, &a)); 228f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatterv(IA, sendcounti, displi, MPIU_INT, ia, n + 1, MPIU_INT, 0, comm)); 229f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatterv(JA, NZ, NZdispl, MPIU_INT, ja, nz, MPIU_INT, 0, comm)); 230f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatterv(sa, NZ, NZdispl, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm)); 231f1f2ae84SBarry Smith 232f1f2ae84SBarry Smith if (pc) { 233f1f2ae84SBarry Smith PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa)); 234f1f2ae84SBarry Smith PetscCall(MatRestoreRowIJ(sA, 0, PETSC_FALSE, PETSC_FALSE, NULL, &IA, &JA, NULL)); 235f1f2ae84SBarry Smith } 236f1f2ae84SBarry Smith 237ad540459SPierre Jolivet for (j = 1; j < n + 1; j++) ia[j] -= ia[0]; 238f1f2ae84SBarry Smith ia[0] = 0; 2390316ec64SBarry Smith PetscCall(PetscLogStagePush(PCMPIStage)); 2403821be0aSBarry Smith PetscCall(MatCreate(comm, &A)); 2413821be0aSBarry Smith if (matproperties[7] > 0) { 2423821be0aSBarry Smith if (!pc) PetscCall(PetscMalloc1(matproperties[7] + 1, &cprefix)); 2433821be0aSBarry Smith PetscCallMPI(MPI_Bcast(cprefix, matproperties[7] + 1, MPI_CHAR, 0, comm)); 2443821be0aSBarry Smith PetscCall(MatSetOptionsPrefix(A, cprefix)); 2453821be0aSBarry Smith PetscCall(PetscFree(cprefix)); 2463821be0aSBarry Smith } 2473821be0aSBarry Smith PetscCall(MatAppendOptionsPrefix(A, "mpi_linear_solver_server_")); 2483821be0aSBarry Smith PetscCall(MatSetSizes(A, m, n, matproperties[0], matproperties[1])); 2493821be0aSBarry Smith PetscCall(MatSetType(A, MATMPIAIJ)); 2503821be0aSBarry Smith PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, a)); 25168a21331SBarry Smith PetscCall(MatSetBlockSize(A, matproperties[2])); 2523821be0aSBarry Smith 25368a21331SBarry Smith if (matproperties[3]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE)); 25468a21331SBarry Smith if (matproperties[4]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[4] == 1 ? PETSC_TRUE : PETSC_FALSE)); 25568a21331SBarry Smith if (matproperties[5]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[5] == 1 ? PETSC_TRUE : PETSC_FALSE)); 25668a21331SBarry Smith if (matproperties[6]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[6] == 1 ? PETSC_TRUE : PETSC_FALSE)); 257f1f2ae84SBarry Smith 258f1f2ae84SBarry Smith PetscCall(PetscFree3(ia, ja, a)); 259f1f2ae84SBarry Smith PetscCall(KSPSetOperators(ksp, A, A)); 260f1f2ae84SBarry Smith if (!ksp->vec_sol) PetscCall(MatCreateVecs(A, &ksp->vec_sol, &ksp->vec_rhs)); 2610316ec64SBarry Smith PetscCall(PetscLogStagePop()); 262f1f2ae84SBarry Smith if (pc) { /* needed for scatterv/gatherv of rhs and solution */ 263f1f2ae84SBarry Smith const PetscInt *range; 264f1f2ae84SBarry Smith 265f1f2ae84SBarry Smith PetscCall(VecGetOwnershipRanges(ksp->vec_sol, &range)); 266f1f2ae84SBarry Smith for (i = 0; i < size; i++) { 267f1f2ae84SBarry Smith km->sendcount[i] = (PetscMPIInt)(range[i + 1] - range[i]); 268f1f2ae84SBarry Smith km->displ[i] = (PetscMPIInt)range[i]; 269f1f2ae84SBarry Smith } 270f1f2ae84SBarry Smith } 271f1f2ae84SBarry Smith PetscCall(MatDestroy(&A)); 272f1f2ae84SBarry Smith PetscCall(KSPSetFromOptions(ksp)); 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274f1f2ae84SBarry Smith } 275f1f2ae84SBarry Smith 276d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIUpdateMatValues(PC pc) 277d71ae5a4SJacob Faibussowitsch { 278f1f2ae84SBarry Smith PC_MPI *km = pc ? (PC_MPI *)pc->data : NULL; 279f1f2ae84SBarry Smith KSP ksp; 280f1f2ae84SBarry Smith Mat sA, A; 281f1f2ae84SBarry Smith MPI_Comm comm = PC_MPI_COMM_WORLD; 282f1f2ae84SBarry Smith PetscScalar *a; 283f1f2ae84SBarry Smith PetscCount nz; 284f1f2ae84SBarry Smith const PetscScalar *sa = NULL; 285dad3da8eSBarry Smith PetscMPIInt size; 28668a21331SBarry Smith PetscInt matproperties[4] = {0, 0, 0, 0}; 287f1f2ae84SBarry Smith 288f1f2ae84SBarry Smith PetscFunctionBegin; 289f1f2ae84SBarry Smith if (pc) { 290f1f2ae84SBarry Smith PetscCall(PCGetOperators(pc, &sA, &sA)); 291f1f2ae84SBarry Smith PetscCall(MatSeqAIJGetArrayRead(sA, &sa)); 292f1f2ae84SBarry Smith } 293f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm)); 2943ba16761SJacob Faibussowitsch if (!ksp) PetscFunctionReturn(PETSC_SUCCESS); 295f1f2ae84SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm)); 296dad3da8eSBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 297dad3da8eSBarry Smith PCMPIMatCounts[size - 1]++; 298f1f2ae84SBarry Smith PetscCall(KSPGetOperators(ksp, NULL, &A)); 299f1f2ae84SBarry Smith PetscCall(MatMPIAIJGetNumberNonzeros(A, &nz)); 300f1f2ae84SBarry Smith PetscCall(PetscMalloc1(nz, &a)); 301f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatterv(sa, pc ? km->NZ : NULL, pc ? km->NZdispl : NULL, MPIU_SCALAR, a, nz, MPIU_SCALAR, 0, comm)); 30268a21331SBarry Smith if (pc) { 30368a21331SBarry Smith PetscBool isset, issymmetric, ishermitian, isspd, isstructurallysymmetric; 30468a21331SBarry Smith 30568a21331SBarry Smith PetscCall(MatSeqAIJRestoreArrayRead(sA, &sa)); 30668a21331SBarry Smith 30768a21331SBarry Smith PetscCall(MatIsSymmetricKnown(sA, &isset, &issymmetric)); 30868a21331SBarry Smith matproperties[0] = !isset ? 0 : (issymmetric ? 1 : 2); 30968a21331SBarry Smith PetscCall(MatIsHermitianKnown(sA, &isset, &ishermitian)); 31068a21331SBarry Smith matproperties[1] = !isset ? 0 : (ishermitian ? 1 : 2); 31168a21331SBarry Smith PetscCall(MatIsSPDKnown(sA, &isset, &isspd)); 31268a21331SBarry Smith matproperties[2] = !isset ? 0 : (isspd ? 1 : 2); 31368a21331SBarry Smith PetscCall(MatIsStructurallySymmetricKnown(sA, &isset, &isstructurallysymmetric)); 31468a21331SBarry Smith matproperties[3] = !isset ? 0 : (isstructurallysymmetric ? 1 : 2); 31568a21331SBarry Smith } 316f1f2ae84SBarry Smith PetscCall(MatUpdateMPIAIJWithArray(A, a)); 317f1f2ae84SBarry Smith PetscCall(PetscFree(a)); 31868a21331SBarry Smith PetscCallMPI(MPI_Bcast(matproperties, 4, MPIU_INT, 0, comm)); 31968a21331SBarry 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 */ 32068a21331SBarry Smith if (matproperties[0]) PetscCall(MatSetOption(A, MAT_SYMMETRIC, matproperties[0] == 1 ? PETSC_TRUE : PETSC_FALSE)); 32168a21331SBarry Smith if (matproperties[1]) PetscCall(MatSetOption(A, MAT_HERMITIAN, matproperties[1] == 1 ? PETSC_TRUE : PETSC_FALSE)); 32268a21331SBarry Smith if (matproperties[2]) PetscCall(MatSetOption(A, MAT_SPD, matproperties[2] == 1 ? PETSC_TRUE : PETSC_FALSE)); 32368a21331SBarry Smith if (matproperties[3]) PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, matproperties[3] == 1 ? PETSC_TRUE : PETSC_FALSE)); 3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 325f1f2ae84SBarry Smith } 326f1f2ae84SBarry Smith 327d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPISolve(PC pc, Vec B, Vec X) 328d71ae5a4SJacob Faibussowitsch { 329f1f2ae84SBarry Smith PC_MPI *km = pc ? (PC_MPI *)pc->data : NULL; 330f1f2ae84SBarry Smith KSP ksp; 331f1f2ae84SBarry Smith MPI_Comm comm = PC_MPI_COMM_WORLD; 332f1f2ae84SBarry Smith const PetscScalar *sb = NULL, *x; 333f1f2ae84SBarry Smith PetscScalar *b, *sx = NULL; 3345316cbedSBarry Smith PetscInt its, n; 3355316cbedSBarry Smith PetscMPIInt size; 336f1f2ae84SBarry Smith 337f1f2ae84SBarry Smith PetscFunctionBegin; 338f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatter(pc ? km->ksps : &ksp, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm)); 3393ba16761SJacob Faibussowitsch if (!ksp) PetscFunctionReturn(PETSC_SUCCESS); 340f1f2ae84SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm)); 341f1f2ae84SBarry Smith 342f1f2ae84SBarry Smith /* TODO: optimize code to not require building counts/displ every time */ 343f1f2ae84SBarry Smith 344f1f2ae84SBarry Smith /* scatterv rhs */ 345dad3da8eSBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 3465316cbedSBarry Smith if (pc) { 3475316cbedSBarry Smith PetscInt N; 3485316cbedSBarry Smith 349dad3da8eSBarry Smith PCMPISolveCounts[size - 1]++; 35068a21331SBarry Smith PetscCall(MatGetSize(pc->pmat, &N, NULL)); 3515316cbedSBarry Smith PCMPISizes[size - 1] += N; 352f1f2ae84SBarry Smith PetscCall(VecGetArrayRead(B, &sb)); 353f1f2ae84SBarry Smith } 354f1f2ae84SBarry Smith PetscCall(VecGetLocalSize(ksp->vec_rhs, &n)); 355f1f2ae84SBarry Smith PetscCall(VecGetArray(ksp->vec_rhs, &b)); 356f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatterv(sb, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, b, n, MPIU_SCALAR, 0, comm)); 357f1f2ae84SBarry Smith PetscCall(VecRestoreArray(ksp->vec_rhs, &b)); 358f1f2ae84SBarry Smith if (pc) PetscCall(VecRestoreArrayRead(B, &sb)); 359f1f2ae84SBarry Smith 3600316ec64SBarry Smith PetscCall(PetscLogStagePush(PCMPIStage)); 361f1f2ae84SBarry Smith PetscCall(KSPSolve(ksp, NULL, NULL)); 3620316ec64SBarry Smith PetscCall(PetscLogStagePop()); 3635316cbedSBarry Smith PetscCall(KSPGetIterationNumber(ksp, &its)); 3645316cbedSBarry Smith PCMPIIterations[size - 1] += its; 365f1f2ae84SBarry Smith 366f1f2ae84SBarry Smith /* gather solution */ 367f1f2ae84SBarry Smith PetscCall(VecGetArrayRead(ksp->vec_sol, &x)); 368f1f2ae84SBarry Smith if (pc) PetscCall(VecGetArray(X, &sx)); 369f1f2ae84SBarry Smith PetscCallMPI(MPI_Gatherv(x, n, MPIU_SCALAR, sx, pc ? km->sendcount : NULL, pc ? km->displ : NULL, MPIU_SCALAR, 0, comm)); 370f1f2ae84SBarry Smith if (pc) PetscCall(VecRestoreArray(X, &sx)); 371f1f2ae84SBarry Smith PetscCall(VecRestoreArrayRead(ksp->vec_sol, &x)); 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373f1f2ae84SBarry Smith } 374f1f2ae84SBarry Smith 375d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMPIDestroy(PC pc) 376d71ae5a4SJacob Faibussowitsch { 377f1f2ae84SBarry Smith PC_MPI *km = pc ? (PC_MPI *)pc->data : NULL; 378f1f2ae84SBarry Smith KSP ksp; 379f1f2ae84SBarry Smith MPI_Comm comm = PC_MPI_COMM_WORLD; 380f1f2ae84SBarry Smith 381f1f2ae84SBarry Smith PetscFunctionBegin; 382f1f2ae84SBarry Smith PetscCallMPI(MPI_Scatter(pc ? km->ksps : NULL, 1, MPI_AINT, &ksp, 1, MPI_AINT, 0, comm)); 3833ba16761SJacob Faibussowitsch if (!ksp) PetscFunctionReturn(PETSC_SUCCESS); 384f1f2ae84SBarry Smith PetscCall(KSPDestroy(&ksp)); 3853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 386f1f2ae84SBarry Smith } 387f1f2ae84SBarry Smith 3883821be0aSBarry Smith PetscBool PCMPIServerActive = PETSC_FALSE; 3893821be0aSBarry Smith 390f1f2ae84SBarry Smith /*@C 3917a99bfcaSBarry Smith PCMPIServerBegin - starts a server that runs on the `rank != 0` MPI processes waiting to process requests for 392f1580f4eSBarry Smith parallel `KSP` solves and management of parallel `KSP` objects. 393f1f2ae84SBarry Smith 3943821be0aSBarry Smith Logically Collective on all MPI processes except rank 0 395f1f2ae84SBarry Smith 396f1580f4eSBarry Smith Options Database Keys: 397f1f2ae84SBarry 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 3983821be0aSBarry Smith - -mpi_linear_solver_server_view - displays information about all the linear systems solved by the MPI linear solver server at the conclusion of the program 399f1f2ae84SBarry Smith 40020f4b53cSBarry Smith Level: developer 40120f4b53cSBarry Smith 402f1580f4eSBarry Smith Note: 403f1f2ae84SBarry Smith This is normally started automatically in `PetscInitialize()` when the option is provided 404f1f2ae84SBarry Smith 4053821be0aSBarry Smith See `PCMPI` for information on using the solver with a `KSP` object 4063821be0aSBarry Smith 407f1f2ae84SBarry Smith Developer Notes: 4083821be0aSBarry Smith When called on MPI rank 0 this sets `PETSC_COMM_WORLD` to `PETSC_COMM_SELF` to allow a main program 409f1580f4eSBarry Smith written with `PETSC_COMM_WORLD` to run correctly on the single rank while all the ranks 410f1580f4eSBarry Smith (that would normally be sharing `PETSC_COMM_WORLD`) to run the solver server. 411f1f2ae84SBarry Smith 412f1580f4eSBarry Smith Can this be integrated into the `PetscDevice` abstraction that is currently being developed? 413f1f2ae84SBarry Smith 4143821be0aSBarry Smith Conceivably `PCREDISTRIBUTE` could be organized in a similar manner to simplify its usage 4153821be0aSBarry Smith 4163821be0aSBarry Smith This could be implemented directly at the `KSP` level instead of using the `PCMPI` wrapper object 4173821be0aSBarry Smith 418baca6076SPierre Jolivet The code could be extended to allow an MPI + OpenMP application to use the linear solver server concept across all shared-memory 4193821be0aSBarry Smith nodes with a single MPI process per node for the user application but multiple MPI processes per node for the linear solver. 4203821be0aSBarry Smith 4217a99bfcaSBarry Smith The concept could also be extended for users's callbacks for `SNES`, `TS`, and `Tao` where the `SNESSolve()` for example, runs on 4223821be0aSBarry Smith all MPI processes but the user callback only runs on one MPI process per node. 4233821be0aSBarry Smith 4243821be0aSBarry Smith PETSc could also be extended with an MPI-less API that provides access to PETSc's solvers without any reference to MPI, essentially remove 4253821be0aSBarry Smith the `MPI_Comm` argument from PETSc calls. 4263821be0aSBarry Smith 4273821be0aSBarry Smith .seealso: [](sec_pcmpi), `PCMPIServerEnd()`, `PCMPI`, `KSPCheckPCMPI()` 428f1f2ae84SBarry Smith @*/ 429d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerBegin(void) 430d71ae5a4SJacob Faibussowitsch { 431f1f2ae84SBarry Smith PetscMPIInt rank; 432f1f2ae84SBarry Smith 433f1f2ae84SBarry Smith PetscFunctionBegin; 4349d3446b2SPierre Jolivet PetscCall(PetscInfo(NULL, "Starting MPI Linear Solver Server\n")); 4355e1a0e3cSBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY)) { 4365e1a0e3cSBarry Smith PetscCall(VecInitializePackage()); 4375e1a0e3cSBarry Smith PetscCall(MatInitializePackage()); 4385e1a0e3cSBarry Smith PetscCall(DMInitializePackage()); 4395e1a0e3cSBarry Smith PetscCall(PCInitializePackage()); 4405e1a0e3cSBarry Smith PetscCall(KSPInitializePackage()); 4415e1a0e3cSBarry Smith PetscCall(SNESInitializePackage()); 4425e1a0e3cSBarry Smith PetscCall(TSInitializePackage()); 4435e1a0e3cSBarry Smith PetscCall(TaoInitializePackage()); 4445e1a0e3cSBarry Smith } 445*956255efSBarry Smith PetscCall(PetscLogStageRegister("PCMPI", &PCMPIStage)); 4465e1a0e3cSBarry Smith 447f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_rank(PC_MPI_COMM_WORLD, &rank)); 448f1f2ae84SBarry Smith if (rank == 0) { 449f1f2ae84SBarry Smith PETSC_COMM_WORLD = PETSC_COMM_SELF; 4503821be0aSBarry Smith PCMPIServerActive = PETSC_TRUE; 4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 452f1f2ae84SBarry Smith } 453f1f2ae84SBarry Smith 454f1f2ae84SBarry Smith while (PETSC_TRUE) { 455f1f2ae84SBarry Smith PCMPICommand request = PCMPI_CREATE; 456f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD)); 457f1f2ae84SBarry Smith switch (request) { 458d71ae5a4SJacob Faibussowitsch case PCMPI_CREATE: 459d71ae5a4SJacob Faibussowitsch PetscCall(PCMPICreate(NULL)); 460d71ae5a4SJacob Faibussowitsch break; 461d71ae5a4SJacob Faibussowitsch case PCMPI_SET_MAT: 462d71ae5a4SJacob Faibussowitsch PetscCall(PCMPISetMat(NULL)); 463d71ae5a4SJacob Faibussowitsch break; 464d71ae5a4SJacob Faibussowitsch case PCMPI_UPDATE_MAT_VALUES: 465d71ae5a4SJacob Faibussowitsch PetscCall(PCMPIUpdateMatValues(NULL)); 466d71ae5a4SJacob Faibussowitsch break; 467f1f2ae84SBarry Smith case PCMPI_VIEW: 468f1f2ae84SBarry Smith // PetscCall(PCMPIView(NULL)); 469f1f2ae84SBarry Smith break; 470d71ae5a4SJacob Faibussowitsch case PCMPI_SOLVE: 471d71ae5a4SJacob Faibussowitsch PetscCall(PCMPISolve(NULL, NULL, NULL)); 472d71ae5a4SJacob Faibussowitsch break; 473d71ae5a4SJacob Faibussowitsch case PCMPI_DESTROY: 474d71ae5a4SJacob Faibussowitsch PetscCall(PCMPIDestroy(NULL)); 475d71ae5a4SJacob Faibussowitsch break; 476f1f2ae84SBarry Smith case PCMPI_EXIT: 477f1f2ae84SBarry Smith PetscCall(PetscFinalize()); 478f1f2ae84SBarry Smith exit(0); /* not sure if this is a good idea, but cannot return because it will run users main program */ 479f1f2ae84SBarry Smith break; 480d71ae5a4SJacob Faibussowitsch default: 481d71ae5a4SJacob Faibussowitsch break; 482f1f2ae84SBarry Smith } 483f1f2ae84SBarry Smith } 4843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 485f1f2ae84SBarry Smith } 486f1f2ae84SBarry Smith 487f1f2ae84SBarry Smith /*@C 488f1f2ae84SBarry Smith PCMPIServerEnd - ends a server that runs on the rank != 0 MPI processes waiting to process requests for 489f1580f4eSBarry Smith parallel KSP solves and management of parallel `KSP` objects. 490f1f2ae84SBarry Smith 49120f4b53cSBarry Smith Logically Collective on all MPI ranks except 0 49220f4b53cSBarry Smith 49320f4b53cSBarry Smith Level: developer 494f1f2ae84SBarry Smith 495f1580f4eSBarry Smith Note: 496f1f2ae84SBarry Smith This is normally ended automatically in `PetscFinalize()` when the option is provided 497f1f2ae84SBarry Smith 4983821be0aSBarry Smith .seealso: [](sec_pcmpi), `PCMPIServerBegin()`, `PCMPI`, `KSPCheckPCMPI()` 499f1f2ae84SBarry Smith @*/ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMPIServerEnd(void) 501d71ae5a4SJacob Faibussowitsch { 502f1f2ae84SBarry Smith PCMPICommand request = PCMPI_EXIT; 503f1f2ae84SBarry Smith 504f1f2ae84SBarry Smith PetscFunctionBegin; 505f1f2ae84SBarry Smith if (PetscGlobalRank == 0) { 506f1f2ae84SBarry Smith PetscViewer viewer = NULL; 507f1f2ae84SBarry Smith PetscViewerFormat format; 508f1f2ae84SBarry Smith 509f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, PC_MPI_COMM_WORLD)); 510f1f2ae84SBarry Smith PETSC_COMM_WORLD = MPI_COMM_WORLD; /* could use PC_MPI_COMM_WORLD */ 511f1f2ae84SBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, NULL, "MPI linear solver server options", NULL); 512f1f2ae84SBarry Smith PetscCall(PetscOptionsViewer("-mpi_linear_solver_server_view", "View information about system solved with the server", "PCMPI", &viewer, &format, NULL)); 513f1f2ae84SBarry Smith PetscOptionsEnd(); 514f1f2ae84SBarry Smith if (viewer) { 515f1f2ae84SBarry Smith PetscBool isascii; 516f1f2ae84SBarry Smith 517f1f2ae84SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 518f1f2ae84SBarry Smith if (isascii) { 519f1f2ae84SBarry Smith PetscMPIInt size; 5205316cbedSBarry Smith PetscMPIInt i; 521f1f2ae84SBarry Smith 522f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 5235316cbedSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "MPI linear solver server statistics:\n")); 5245316cbedSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Ranks KSPSolve()s Mats KSPs Avg. Size Avg. Its\n")); 5255316cbedSBarry Smith if (PCMPIKSPCountsSeq) { 5265316cbedSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Sequential %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", PCMPISolveCountsSeq, PCMPIKSPCountsSeq, PCMPISizesSeq / PCMPISolveCountsSeq, PCMPIIterationsSeq / PCMPISolveCountsSeq)); 527f1f2ae84SBarry Smith } 5285316cbedSBarry Smith for (i = 0; i < size; i++) { 5295316cbedSBarry Smith if (PCMPIKSPCounts[i]) { 5305316cbedSBarry 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])); 5315316cbedSBarry Smith } 5325316cbedSBarry Smith } 533f1f2ae84SBarry Smith } 534f1f2ae84SBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 535f1f2ae84SBarry Smith } 536f1f2ae84SBarry Smith } 537f1f2ae84SBarry Smith PetscCall(PCMPICommsDestroy()); 5383821be0aSBarry Smith PCMPIServerActive = PETSC_FALSE; 5393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 540f1f2ae84SBarry Smith } 541f1f2ae84SBarry Smith 542f1f2ae84SBarry Smith /* 543f1f2ae84SBarry Smith This version is used in the trivial case when the MPI parallel solver server is running on just the original MPI rank 0 544f1f2ae84SBarry Smith because, for example, the problem is small. This version is more efficient because it does not require copying any data 545f1f2ae84SBarry Smith */ 546d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Seq(PC pc) 547d71ae5a4SJacob Faibussowitsch { 548f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 549f1f2ae84SBarry Smith Mat sA; 550f1f2ae84SBarry Smith const char *prefix; 551f9818f3cSJose E. Roman char *found = NULL, *cprefix; 552f1f2ae84SBarry Smith 553f1f2ae84SBarry Smith PetscFunctionBegin; 554f1f2ae84SBarry Smith PetscCall(PCGetOperators(pc, NULL, &sA)); 555f1f2ae84SBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 556f1f2ae84SBarry Smith PetscCall(KSPCreate(PETSC_COMM_SELF, &km->ksps[0])); 5573821be0aSBarry Smith PetscCall(KSPSetNestLevel(km->ksps[0], 1)); 5583821be0aSBarry Smith PetscCall(PetscObjectSetTabLevel((PetscObject)km->ksps[0], 1)); 5593821be0aSBarry Smith 5603821be0aSBarry Smith /* Created KSP gets prefix of PC minus the mpi_linear_solver_server_ portion */ 5613821be0aSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 5623821be0aSBarry Smith PetscCheck(prefix, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing required prefix"); 5633821be0aSBarry Smith PetscCall(PetscStrallocpy(prefix, &cprefix)); 5643821be0aSBarry Smith PetscCall(PetscStrstr(cprefix, "mpi_linear_solver_server_", &found)); 5653821be0aSBarry Smith PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI missing mpi_linear_solver_server_ portion of prefix"); 5663821be0aSBarry Smith *found = 0; 5673821be0aSBarry Smith PetscCall(KSPSetOptionsPrefix(km->ksps[0], cprefix)); 5683821be0aSBarry Smith PetscCall(PetscFree(cprefix)); 5693821be0aSBarry Smith 570f1f2ae84SBarry Smith PetscCall(KSPSetOperators(km->ksps[0], sA, sA)); 571f1f2ae84SBarry Smith PetscCall(KSPSetFromOptions(km->ksps[0])); 572f1f2ae84SBarry Smith PetscCall(KSPSetUp(km->ksps[0])); 5733ba16761SJacob Faibussowitsch PetscCall(PetscInfo((PetscObject)pc, "MPI parallel linear solver system is being solved directly on rank 0 due to its small size\n")); 574f1f2ae84SBarry Smith PCMPIKSPCountsSeq++; 5753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 576f1f2ae84SBarry Smith } 577f1f2ae84SBarry Smith 578d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Seq(PC pc, Vec b, Vec x) 579d71ae5a4SJacob Faibussowitsch { 580f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 5815316cbedSBarry Smith PetscInt its, n; 5825316cbedSBarry Smith Mat A; 583f1f2ae84SBarry Smith 584f1f2ae84SBarry Smith PetscFunctionBegin; 585f1f2ae84SBarry Smith PetscCall(KSPSolve(km->ksps[0], b, x)); 5865316cbedSBarry Smith PetscCall(KSPGetIterationNumber(km->ksps[0], &its)); 587f1f2ae84SBarry Smith PCMPISolveCountsSeq++; 5885316cbedSBarry Smith PCMPIIterationsSeq += its; 5895316cbedSBarry Smith PetscCall(KSPGetOperators(km->ksps[0], NULL, &A)); 5905316cbedSBarry Smith PetscCall(MatGetSize(A, &n, NULL)); 5915316cbedSBarry Smith PCMPISizesSeq += n; 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593f1f2ae84SBarry Smith } 594f1f2ae84SBarry Smith 595d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Seq(PC pc, PetscViewer viewer) 596d71ae5a4SJacob Faibussowitsch { 597f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 598f1f2ae84SBarry Smith 599f1f2ae84SBarry Smith PetscFunctionBegin; 600f1f2ae84SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Running MPI linear solver server directly on rank 0 due to its small size\n")); 601f1f2ae84SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank)); 6025316cbedSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n")); 6033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 604f1f2ae84SBarry Smith } 605f1f2ae84SBarry Smith 606d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Seq(PC pc) 607d71ae5a4SJacob Faibussowitsch { 608f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 609f1f2ae84SBarry Smith 610f1f2ae84SBarry Smith PetscFunctionBegin; 611f1f2ae84SBarry Smith PetscCall(KSPDestroy(&km->ksps[0])); 612f1f2ae84SBarry Smith PetscCall(PetscFree(pc->data)); 6133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 614f1f2ae84SBarry Smith } 615f1f2ae84SBarry Smith 616f1f2ae84SBarry Smith /* 617f1f2ae84SBarry Smith PCSetUp_MPI - Trigger the creation of the MPI parallel PC and copy parts of the matrix and 618f1f2ae84SBarry Smith right hand side to the parallel PC 619f1f2ae84SBarry Smith */ 620d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_MPI(PC pc) 621d71ae5a4SJacob Faibussowitsch { 622f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 623f1f2ae84SBarry Smith PetscMPIInt rank, size; 624f1f2ae84SBarry Smith PCMPICommand request; 625f1f2ae84SBarry Smith PetscBool newmatrix = PETSC_FALSE; 626f1f2ae84SBarry Smith 627f1f2ae84SBarry Smith PetscFunctionBegin; 628f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank)); 629f1f2ae84SBarry 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?"); 630f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size)); 631f1f2ae84SBarry Smith 632f1f2ae84SBarry Smith if (!pc->setupcalled) { 633f1f2ae84SBarry Smith if (!km->alwaysuseserver) { 634f1f2ae84SBarry Smith PetscInt n; 635f1f2ae84SBarry Smith Mat sA; 636f1f2ae84SBarry Smith /* short circuit for small systems */ 637f1f2ae84SBarry Smith PetscCall(PCGetOperators(pc, &sA, &sA)); 638f1f2ae84SBarry Smith PetscCall(MatGetSize(sA, &n, NULL)); 639f1f2ae84SBarry Smith if (n < 2 * km->mincntperrank - 1 || size == 1) { 640f1f2ae84SBarry Smith pc->ops->setup = NULL; 641f1f2ae84SBarry Smith pc->ops->apply = PCApply_Seq; 642f1f2ae84SBarry Smith pc->ops->destroy = PCDestroy_Seq; 643f1f2ae84SBarry Smith pc->ops->view = PCView_Seq; 644f1f2ae84SBarry Smith PetscCall(PCSetUp_Seq(pc)); 6453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 646f1f2ae84SBarry Smith } 647f1f2ae84SBarry Smith } 648f1f2ae84SBarry Smith 649f1f2ae84SBarry Smith request = PCMPI_CREATE; 650f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD)); 651f1f2ae84SBarry Smith PetscCall(PCMPICreate(pc)); 652f1f2ae84SBarry Smith newmatrix = PETSC_TRUE; 6539371c9d4SSatish Balay } 6549371c9d4SSatish Balay if (pc->flag == DIFFERENT_NONZERO_PATTERN) newmatrix = PETSC_TRUE; 655f1f2ae84SBarry Smith 656f1f2ae84SBarry Smith if (newmatrix) { 6573ba16761SJacob Faibussowitsch PetscCall(PetscInfo((PetscObject)pc, "New matrix or matrix has changed nonzero structure\n")); 658f1f2ae84SBarry Smith request = PCMPI_SET_MAT; 659f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD)); 660f1f2ae84SBarry Smith PetscCall(PCMPISetMat(pc)); 661f1f2ae84SBarry Smith } else { 662bbea24aaSStefano Zampini PetscCall(PetscInfo((PetscObject)pc, "Matrix has only changed nonzero values\n")); 663f1f2ae84SBarry Smith request = PCMPI_UPDATE_MAT_VALUES; 664f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD)); 665f1f2ae84SBarry Smith PetscCall(PCMPIUpdateMatValues(pc)); 666f1f2ae84SBarry Smith } 6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 668f1f2ae84SBarry Smith } 669f1f2ae84SBarry Smith 670d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MPI(PC pc, Vec b, Vec x) 671d71ae5a4SJacob Faibussowitsch { 672f1f2ae84SBarry Smith PCMPICommand request = PCMPI_SOLVE; 673f1f2ae84SBarry Smith 674f1f2ae84SBarry Smith PetscFunctionBegin; 675f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD)); 676f1f2ae84SBarry Smith PetscCall(PCMPISolve(pc, b, x)); 6773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 678f1f2ae84SBarry Smith } 679f1f2ae84SBarry Smith 68066976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_MPI(PC pc) 681d71ae5a4SJacob Faibussowitsch { 682f1f2ae84SBarry Smith PCMPICommand request = PCMPI_DESTROY; 683f1f2ae84SBarry Smith 684f1f2ae84SBarry Smith PetscFunctionBegin; 685f1f2ae84SBarry Smith PetscCallMPI(MPI_Bcast(&request, 1, MPIU_ENUM, 0, MPI_COMM_WORLD)); 686f1f2ae84SBarry Smith PetscCall(PCMPIDestroy(pc)); 687f1f2ae84SBarry Smith PetscCall(PetscFree(pc->data)); 6883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 689f1f2ae84SBarry Smith } 690f1f2ae84SBarry Smith 691f1f2ae84SBarry Smith /* 692f1f2ae84SBarry Smith PCView_MPI - Cannot call view on the MPI parallel KSP because other ranks do not have access to the viewer 693f1f2ae84SBarry Smith */ 69466976f2fSJacob Faibussowitsch static PetscErrorCode PCView_MPI(PC pc, PetscViewer viewer) 695d71ae5a4SJacob Faibussowitsch { 696f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 697f1f2ae84SBarry Smith MPI_Comm comm; 698f1f2ae84SBarry Smith PetscMPIInt size; 699f1f2ae84SBarry Smith 700f1f2ae84SBarry Smith PetscFunctionBegin; 701f1f2ae84SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)km->ksps[0], &comm)); 702f1f2ae84SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 703f1f2ae84SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Size of MPI communicator used for MPI parallel KSP solve %d\n", size)); 704f1f2ae84SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Desired minimum number of nonzeros per rank for MPI parallel solve %d\n", (int)km->mincntperrank)); 70568a21331SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "*** Use -mpi_linear_solver_server_view to statistics on all the solves ***\n")); 7063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 707f1f2ae84SBarry Smith } 708f1f2ae84SBarry Smith 70966976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_MPI(PC pc, PetscOptionItems *PetscOptionsObject) 710d71ae5a4SJacob Faibussowitsch { 711f1f2ae84SBarry Smith PC_MPI *km = (PC_MPI *)pc->data; 712f1f2ae84SBarry Smith 713f1f2ae84SBarry Smith PetscFunctionBegin; 714f1f2ae84SBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "MPI linear solver server options"); 7153821be0aSBarry Smith PetscCall(PetscOptionsInt("-minimum_count_per_rank", "Desired minimum number of nonzeros per rank", "None", km->mincntperrank, &km->mincntperrank, NULL)); 7163821be0aSBarry Smith PetscCall(PetscOptionsBool("-always_use_server", "Use the server even if only one rank is used for the solve (for debugging)", "None", km->alwaysuseserver, &km->alwaysuseserver, NULL)); 717f1f2ae84SBarry Smith PetscOptionsHeadEnd(); 7183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 719f1f2ae84SBarry Smith } 720f1f2ae84SBarry Smith 721f1f2ae84SBarry Smith /*MC 722f1580f4eSBarry Smith PCMPI - Calls an MPI parallel `KSP` to solve a linear system from user code running on one process 723f1f2ae84SBarry Smith 7243821be0aSBarry Smith Options Database Keys for the Server: 725f1f2ae84SBarry 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 7263821be0aSBarry Smith - -mpi_linear_solver_server_view - displays information about all the linear systems solved by the MPI linear solver server 727f1f2ae84SBarry Smith 7283821be0aSBarry Smith Options Database Keys for a specific `KSP` object 7293821be0aSBarry Smith + -[any_ksp_prefix]_mpi_linear_solver_server_minimum_count_per_rank - sets the minimum size of the linear system per MPI rank that the solver will strive for 7303821be0aSBarry Smith - -[any_ksp_prefix]_mpi_linear_solver_server_always_use_server - use the server solver code even if the particular system is only solved on the process (for debugging and testing purposes) 7313821be0aSBarry Smith 7323821be0aSBarry Smith Level: developer 73320f4b53cSBarry Smith 734f1f2ae84SBarry Smith Notes: 73546bbbc36SPierre Jolivet The options database prefix for the actual solver is any prefix provided before use to the original `KSP` with `KSPSetOptionsPrefix()`, mostly commonly no prefix is used. 736f1f2ae84SBarry Smith 737f1f2ae84SBarry Smith It can be particularly useful for user OpenMP code or potentially user GPU code. 738f1f2ae84SBarry Smith 7393821be0aSBarry Smith When the program is running with a single MPI process then it directly uses the provided matrix and right hand side 7403821be0aSBarry Smith and does not need to distribute the matrix and vector to the various MPI processes; thus it incurs no extra overhead over just using the `KSP` directly. 741f1f2ae84SBarry Smith 7423821be0aSBarry Smith The solver options for actual solving `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 7430316ec64SBarry Smith because they are not the actual solver objects. 7440316ec64SBarry Smith 7450316ec64SBarry 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 74668a21331SBarry 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. 7470316ec64SBarry Smith 7483821be0aSBarry Smith Developer Note: 7493821be0aSBarry Smith This `PCType` is never directly selected by the user, it is set when the option `-mpi_linear_solver_server` is used and the `PC` is at the outer most nesting of 7503821be0aSBarry Smith a `KSP`. The outer most `KSP` object is automatically set to `KSPPREONLY` and thus is not directly visible to the user. 7513821be0aSBarry Smith 7523821be0aSBarry Smith .seealso: [](sec_pcmpi), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `PC`, `PCMPIServerBegin()`, `PCMPIServerEnd()`, `KSPCheckPCMPI()` 753f1f2ae84SBarry Smith M*/ 754d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MPI(PC pc) 755d71ae5a4SJacob Faibussowitsch { 756f1f2ae84SBarry Smith PC_MPI *km; 757f9818f3cSJose E. Roman char *found = NULL; 758f1f2ae84SBarry Smith 759f1f2ae84SBarry Smith PetscFunctionBegin; 7603821be0aSBarry Smith PetscCall(PetscStrstr(((PetscObject)pc)->prefix, "mpi_linear_solver_server_", &found)); 7613821be0aSBarry Smith PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCMPI object prefix does not have mpi_linear_solver_server_"); 7623821be0aSBarry Smith 7633821be0aSBarry Smith /* material from PCSetType() */ 7643821be0aSBarry Smith PetscTryTypeMethod(pc, destroy); 7653821be0aSBarry Smith pc->ops->destroy = NULL; 7663821be0aSBarry Smith pc->data = NULL; 7673821be0aSBarry Smith 7683821be0aSBarry Smith PetscCall(PetscFunctionListDestroy(&((PetscObject)pc)->qlist)); 7693821be0aSBarry Smith PetscCall(PetscMemzero(pc->ops, sizeof(struct _PCOps))); 7703821be0aSBarry Smith pc->modifysubmatrices = NULL; 7713821be0aSBarry Smith pc->modifysubmatricesP = NULL; 7723821be0aSBarry Smith pc->setupcalled = 0; 7733821be0aSBarry Smith 7744dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&km)); 775f1f2ae84SBarry Smith pc->data = (void *)km; 776f1f2ae84SBarry Smith 777f1f2ae84SBarry Smith km->mincntperrank = 10000; 778f1f2ae84SBarry Smith 779f1f2ae84SBarry Smith pc->ops->setup = PCSetUp_MPI; 780f1f2ae84SBarry Smith pc->ops->apply = PCApply_MPI; 781f1f2ae84SBarry Smith pc->ops->destroy = PCDestroy_MPI; 782f1f2ae84SBarry Smith pc->ops->view = PCView_MPI; 783f1f2ae84SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MPI; 7843821be0aSBarry Smith PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCMPI)); 7853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 786f1f2ae84SBarry Smith } 787