1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Test PC redistribute on matrix with load imbalance. \n\ 3c4762a1bSJed Brown Modified from src/ksp/ksp/tutorials/ex2.c.\n\ 4c4762a1bSJed Brown Input parameters include:\n\ 5c4762a1bSJed Brown -random_exact_sol : use a random exact solution vector\n\ 6c4762a1bSJed Brown -view_exact_sol : write exact solution vector to stdout\n\ 710582eeaSJed Brown -n <mesh_y> : number of mesh points\n\n"; 8c4762a1bSJed Brown /* 9c4762a1bSJed Brown Example: 10c4762a1bSJed Brown mpiexec -n 8 ./ex3 -n 10000 -ksp_type cg -pc_type bjacobi -sub_pc_type icc -ksp_rtol 1.e-8 -log_view 11c4762a1bSJed Brown mpiexec -n 8 ./ex3 -n 10000 -ksp_type preonly -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc -redistribute_ksp_rtol 1.e-8 -log_view 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include <petscksp.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown int main(int argc,char **args) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown Vec x,b,u; /* approx solution, RHS, exact solution */ 19c4762a1bSJed Brown Mat A; /* linear system matrix */ 20c4762a1bSJed Brown KSP ksp; /* linear solver context */ 21c4762a1bSJed Brown PetscRandom rctx; /* random number generator context */ 22c4762a1bSJed Brown PetscReal norm; /* norm of solution error */ 23c4762a1bSJed Brown PetscInt i,j,Ii,J,Istart,Iend,m,n = 7,its,nloc,matdistribute=0; 24c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 25c4762a1bSJed Brown PetscScalar v; 26c4762a1bSJed Brown PetscMPIInt rank,size; 27c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 28c4762a1bSJed Brown PetscLogStage stage; 29c4762a1bSJed Brown #endif 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 347827d75bSBarry Smith PetscCheck(size > 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least 2 MPI processes!"); 35c4762a1bSJed Brown 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-matdistribute",&matdistribute,NULL)); 38c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 39c4762a1bSJed Brown Compute the matrix and right-hand-side vector that define 40c4762a1bSJed Brown the linear system, Ax = b. 41c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 42c4762a1bSJed Brown switch(matdistribute) { 43c4762a1bSJed Brown case 1: /* very imbalanced process load for matrix A */ 44c4762a1bSJed Brown m = (1+size)*size; 45c4762a1bSJed Brown nloc = (rank+1)*n; 46c4762a1bSJed Brown if (rank == size-1) { /* proc[size-1] stores all remaining rows */ 47c4762a1bSJed Brown nloc = m*n; 48c4762a1bSJed Brown for (i=0; i<size-1; i++) { 49c4762a1bSJed Brown nloc -= (i+1)*n; 50c4762a1bSJed Brown } 51c4762a1bSJed Brown } 52c4762a1bSJed Brown break; 53c4762a1bSJed Brown default: /* proc[0] and proc[1] load much smaller row blocks, the rest processes have same loads */ 54c4762a1bSJed Brown if (rank == 0 || rank == 1) { 55c4762a1bSJed Brown nloc = n; 56c4762a1bSJed Brown } else { 57c4762a1bSJed Brown nloc = 10*n; /* 10x larger load */ 58c4762a1bSJed Brown } 59c4762a1bSJed Brown m = 2 + (size-2)*10; 60c4762a1bSJed Brown break; 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,nloc,nloc,PETSC_DECIDE,PETSC_DECIDE)); 649566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 659566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 669566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,5,NULL)); 679566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 70c4762a1bSJed Brown nloc = Iend-Istart; 719566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] A Istart,Iend: %D %D; nloc %D\n",rank,Istart,Iend,nloc)); 729566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Assembly", &stage)); 759566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 76c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 77c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 789566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 799566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 809566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 819566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 829566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 83c4762a1bSJed Brown } 849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 869566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */ 899566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Create parallel vectors. */ 929566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&u)); 939566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u,nloc,PETSC_DECIDE)); 949566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u)); 959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&b)); 969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&x)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* Set exact solution; then compute right-hand-side vector. */ 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-random_exact_sol",&flg,NULL)); 100c4762a1bSJed Brown if (flg) { 1019566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 1029566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 1039566063dSJacob Faibussowitsch PetscCall(VecSetRandom(u,rctx)); 1049566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 105c4762a1bSJed Brown } else { 1069566063dSJacob Faibussowitsch PetscCall(VecSet(u,1.0)); 107c4762a1bSJed Brown } 1089566063dSJacob Faibussowitsch PetscCall(MatMult(A,u,b)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* View the exact solution vector if desired */ 111c4762a1bSJed Brown flg = PETSC_FALSE; 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL)); 1139566063dSJacob Faibussowitsch if (flg) PetscCall(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Create the linear solver and set various options 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1189566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp)); 1199566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp,A,A)); 120*d0609cedSBarry Smith PetscCall(KSPSetTolerances(ksp,1.e-2/((m+1)*(n+1)),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 1219566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Solve the linear system 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1269566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp,b,x)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Check solution and clean up 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1319566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,-1.0,u)); 1329566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm)); 1339566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp,&its)); 1349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %D\n",(double)norm,its)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* Free work space. */ 1379566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 1389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); PetscCall(VecDestroy(&x)); 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); PetscCall(MatDestroy(&A)); 1409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 141b122ec5aSJacob Faibussowitsch return 0; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown /*TEST 145c4762a1bSJed Brown 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown nsize: 8 148c4762a1bSJed Brown args: -n 100 -ksp_type cg -pc_type bjacobi -sub_pc_type icc -ksp_rtol 1.e-8 149c4762a1bSJed Brown 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: 2 152c4762a1bSJed Brown nsize: 8 153c4762a1bSJed Brown args: -n 100 -ksp_type preonly -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc -redistribute_ksp_rtol 1.e-8 154c4762a1bSJed Brown 155c4762a1bSJed Brown TEST*/ 156