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 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Create and initialize vectors */ 249566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b)); 259566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&ustar)); 269566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u)); 279566063dSJacob Faibussowitsch PetscCall(VecSet(ustar,1.0)); 289566063dSJacob Faibussowitsch PetscCall(VecSet(u,0.0)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Create and assemble matrix */ 319566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&mat)); 329566063dSJacob Faibussowitsch PetscCall(MatSetType(mat,MATSEQAIJ)); 339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat,n,n,n,n)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 359566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat,3,NULL)); 369566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(mat,1,3,NULL)); 379566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(mat,1,3,NULL)); 389566063dSJacob Faibussowitsch PetscCall(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; 429566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown i = n - 1; col[0] = n - 2; col[1] = n - 1; 459566063dSJacob Faibussowitsch PetscCall(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; 479566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES)); 489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Compute right-hand-side vector */ 529566063dSJacob Faibussowitsch PetscCall(MatMult(mat,ustar,b)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Create PC context and set up data structures */ 559566063dSJacob Faibussowitsch PetscCall(PCCreate(PETSC_COMM_WORLD,&pc)); 569566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 579566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); 589566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc,mat,mat)); 599566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Create KSP context and set up data structures */ 629566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp)); 639566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPRICHARDSON)); 649566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 659566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc,mat,mat)); 669566063dSJacob Faibussowitsch PetscCall(KSPSetPC(ksp,pc)); 679566063dSJacob Faibussowitsch PetscCall(KSPSetUp(ksp)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Solve the problem */ 709566063dSJacob Faibussowitsch PetscCall(KSPGetType(ksp,&kspname)); 719566063dSJacob Faibussowitsch PetscCall(PCGetType(pc,&pcname)); 729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname)); 739566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp,b,u)); 749566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp,&its)); 75*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of iterations %" PetscInt_FMT "\n",its)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Free data structures */ 789566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ustar)); 819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 839566063dSJacob Faibussowitsch PetscCall(PCDestroy(&pc)); 849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 85b122ec5aSJacob 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