1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Demonstrates the use of fast Richardson for SOR. And\n\ 3c4762a1bSJed Brown also tests the MatSOR() routines. Input parameters are:\n\ 4c4762a1bSJed Brown -n <n> : problem dimension\n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petscksp.h> 7c4762a1bSJed Brown #include <petscpc.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown int main(int argc,char **args) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown Mat mat; /* matrix */ 12c4762a1bSJed Brown Vec b,ustar,u; /* vectors (RHS, exact solution, approx solution) */ 13c4762a1bSJed Brown PC pc; /* PC context */ 14c4762a1bSJed Brown KSP ksp; /* KSP context */ 15c4762a1bSJed Brown PetscInt n = 10,i,its,col[3]; 16c4762a1bSJed Brown PetscScalar value[3]; 17c4762a1bSJed Brown KSPType kspname; 18c4762a1bSJed Brown PCType pcname; 19c4762a1bSJed Brown 20*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Create and initialize vectors */ 245f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&b)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&ustar)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&u)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(ustar,1.0)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u,0.0)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Create and assemble matrix */ 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&mat)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(mat,MATSEQAIJ)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(mat,n,n,n,n)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(mat)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(mat,3,NULL)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(mat,1,3,NULL)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSBAIJSetPreallocation(mat,1,3,NULL)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSELLSetPreallocation(mat,3,NULL)); 39c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 40c4762a1bSJed Brown for (i=1; i<n-1; i++) { 41c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown i = n - 1; col[0] = n - 2; col[1] = n - 1; 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES)); 46c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Compute right-hand-side vector */ 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,ustar,b)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Create PC context and set up data structures */ 555f80ce2aSJacob Faibussowitsch CHKERRQ(PCCreate(PETSC_COMM_WORLD,&pc)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCNONE)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetFromOptions(pc)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOperators(pc,mat,mat)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetUp(pc)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Create KSP context and set up data structures */ 625f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&ksp)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(ksp,KSPRICHARDSON)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(ksp)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOperators(pc,mat,mat)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetPC(ksp,pc)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetUp(ksp)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Solve the problem */ 705f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetType(ksp,&kspname)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetType(pc,&pcname)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(ksp,b,u)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetIterationNumber(ksp,&its)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of iterations %D\n",its)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Free data structures */ 785f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&ksp)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ustar)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PCDestroy(&pc)); 84*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 85*b122ec5aSJacob Faibussowitsch return 0; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /*TEST 89c4762a1bSJed Brown 90c4762a1bSJed Brown testset: 91c4762a1bSJed Brown args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric 92c4762a1bSJed Brown output_file: output/ex3_1.out 93c4762a1bSJed Brown test: 94c4762a1bSJed Brown suffix: sor_aij 95c4762a1bSJed Brown test: 96c4762a1bSJed Brown suffix: sor_seqbaij 97c4762a1bSJed Brown args: -mat_type seqbaij 98c4762a1bSJed Brown test: 99c4762a1bSJed Brown suffix: sor_seqsbaij 100c4762a1bSJed Brown args: -mat_type seqbaij 101c4762a1bSJed Brown test: 102c4762a1bSJed Brown suffix: sor_seqsell 103c4762a1bSJed Brown args: -mat_type seqsell 104c4762a1bSJed Brown 105c4762a1bSJed Brown TEST*/ 106