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 PetscErrorCode ierr; 25c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 26c4762a1bSJed Brown PetscScalar v; 27c4762a1bSJed Brown PetscMPIInt rank,size; 28c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 29c4762a1bSJed Brown PetscLogStage stage; 30c4762a1bSJed Brown #endif 31c4762a1bSJed Brown 32*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 352c71b3e2SJacob Faibussowitsch PetscCheckFalse(size < 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least 2 MPI processes!"); 36c4762a1bSJed Brown 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-matdistribute",&matdistribute,NULL)); 39c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 40c4762a1bSJed Brown Compute the matrix and right-hand-side vector that define 41c4762a1bSJed Brown the linear system, Ax = b. 42c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 43c4762a1bSJed Brown switch(matdistribute) { 44c4762a1bSJed Brown case 1: /* very imbalanced process load for matrix A */ 45c4762a1bSJed Brown m = (1+size)*size; 46c4762a1bSJed Brown nloc = (rank+1)*n; 47c4762a1bSJed Brown if (rank == size-1) { /* proc[size-1] stores all remaining rows */ 48c4762a1bSJed Brown nloc = m*n; 49c4762a1bSJed Brown for (i=0; i<size-1; i++) { 50c4762a1bSJed Brown nloc -= (i+1)*n; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53c4762a1bSJed Brown break; 54c4762a1bSJed Brown default: /* proc[0] and proc[1] load much smaller row blocks, the rest processes have same loads */ 55c4762a1bSJed Brown if (rank == 0 || rank == 1) { 56c4762a1bSJed Brown nloc = n; 57c4762a1bSJed Brown } else { 58c4762a1bSJed Brown nloc = 10*n; /* 10x larger load */ 59c4762a1bSJed Brown } 60c4762a1bSJed Brown m = 2 + (size-2)*10; 61c4762a1bSJed Brown break; 62c4762a1bSJed Brown } 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,nloc,nloc,PETSC_DECIDE,PETSC_DECIDE)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 69c4762a1bSJed Brown 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend)); 71c4762a1bSJed Brown nloc = Iend-Istart; 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] A Istart,Iend: %D %D; nloc %D\n",rank,Istart,Iend,nloc)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 74c4762a1bSJed Brown 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Assembly", &stage)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 77c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 78c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 795f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 805f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 815f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 825f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 835f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 84c4762a1bSJed Brown } 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */ 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Create parallel vectors. */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&u)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(u,nloc,PETSC_DECIDE)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(u)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&b)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b,&x)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Set exact solution; then compute right-hand-side vector. */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-random_exact_sol",&flg,NULL)); 101c4762a1bSJed Brown if (flg) { 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(u,rctx)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 106c4762a1bSJed Brown } else { 1075f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u,1.0)); 108c4762a1bSJed Brown } 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,u,b)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* View the exact solution vector if desired */ 112c4762a1bSJed Brown flg = PETSC_FALSE; 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL)); 1145f80ce2aSJacob Faibussowitsch if (flg) CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Create the linear solver and set various options 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1195f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&ksp)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ksp,A,A)); 121b8582c7cSSatish Balay ierr = KSPSetTolerances(ksp,1.e-2/((m+1)*(n+1)),PETSC_DEFAULT,PETSC_DEFAULT, 122c4762a1bSJed Brown PETSC_DEFAULT);CHKERRQ(ierr); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(ksp)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126c4762a1bSJed Brown Solve the linear system 127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1285f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(ksp,b,x)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131c4762a1bSJed Brown Check solution and clean up 132c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1335f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,-1.0,u)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetIterationNumber(ksp,&its)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %D\n",(double)norm,its)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Free work space. */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&ksp)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); CHKERRQ(VecDestroy(&x)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); CHKERRQ(MatDestroy(&A)); 142*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 143*b122ec5aSJacob Faibussowitsch return 0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /*TEST 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown nsize: 8 150c4762a1bSJed Brown args: -n 100 -ksp_type cg -pc_type bjacobi -sub_pc_type icc -ksp_rtol 1.e-8 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: 2 154c4762a1bSJed Brown nsize: 8 155c4762a1bSJed 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 156c4762a1bSJed Brown 157c4762a1bSJed Brown TEST*/ 158