1 2 static char help[] = "Demonstrates the use of fast Richardson for SOR. And\n\ 3 also tests the MatSOR() routines. Input parameters are:\n\ 4 -n <n> : problem dimension\n\n"; 5 6 #include <petscksp.h> 7 #include <petscpc.h> 8 9 int main(int argc,char **args) 10 { 11 Mat mat; /* matrix */ 12 Vec b,ustar,u; /* vectors (RHS, exact solution, approx solution) */ 13 PC pc; /* PC context */ 14 KSP ksp; /* KSP context */ 15 PetscInt n = 10,i,its,col[3]; 16 PetscScalar value[3]; 17 KSPType kspname; 18 PCType pcname; 19 20 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 21 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 22 23 /* Create and initialize vectors */ 24 PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b)); 25 PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&ustar)); 26 PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u)); 27 PetscCall(VecSet(ustar,1.0)); 28 PetscCall(VecSet(u,0.0)); 29 30 /* Create and assemble matrix */ 31 PetscCall(MatCreate(PETSC_COMM_SELF,&mat)); 32 PetscCall(MatSetType(mat,MATSEQAIJ)); 33 PetscCall(MatSetSizes(mat,n,n,n,n)); 34 PetscCall(MatSetFromOptions(mat)); 35 PetscCall(MatSeqAIJSetPreallocation(mat,3,NULL)); 36 PetscCall(MatSeqBAIJSetPreallocation(mat,1,3,NULL)); 37 PetscCall(MatSeqSBAIJSetPreallocation(mat,1,3,NULL)); 38 PetscCall(MatSeqSELLSetPreallocation(mat,3,NULL)); 39 value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 40 for (i=1; i<n-1; i++) { 41 col[0] = i-1; col[1] = i; col[2] = i+1; 42 PetscCall(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES)); 43 } 44 i = n - 1; col[0] = n - 2; col[1] = n - 1; 45 PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES)); 46 i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 47 PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES)); 48 PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 49 PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 50 51 /* Compute right-hand-side vector */ 52 PetscCall(MatMult(mat,ustar,b)); 53 54 /* Create PC context and set up data structures */ 55 PetscCall(PCCreate(PETSC_COMM_WORLD,&pc)); 56 PetscCall(PCSetType(pc,PCNONE)); 57 PetscCall(PCSetFromOptions(pc)); 58 PetscCall(PCSetOperators(pc,mat,mat)); 59 PetscCall(PCSetUp(pc)); 60 61 /* Create KSP context and set up data structures */ 62 PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp)); 63 PetscCall(KSPSetType(ksp,KSPRICHARDSON)); 64 PetscCall(KSPSetFromOptions(ksp)); 65 PetscCall(PCSetOperators(pc,mat,mat)); 66 PetscCall(KSPSetPC(ksp,pc)); 67 PetscCall(KSPSetUp(ksp)); 68 69 /* Solve the problem */ 70 PetscCall(KSPGetType(ksp,&kspname)); 71 PetscCall(PCGetType(pc,&pcname)); 72 PetscCall(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname)); 73 PetscCall(KSPSolve(ksp,b,u)); 74 PetscCall(KSPGetIterationNumber(ksp,&its)); 75 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of iterations %D\n",its)); 76 77 /* Free data structures */ 78 PetscCall(KSPDestroy(&ksp)); 79 PetscCall(VecDestroy(&u)); 80 PetscCall(VecDestroy(&ustar)); 81 PetscCall(VecDestroy(&b)); 82 PetscCall(MatDestroy(&mat)); 83 PetscCall(PCDestroy(&pc)); 84 PetscCall(PetscFinalize()); 85 return 0; 86 } 87 88 /*TEST 89 90 testset: 91 args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric 92 output_file: output/ex3_1.out 93 test: 94 suffix: sor_aij 95 test: 96 suffix: sor_seqbaij 97 args: -mat_type seqbaij 98 test: 99 suffix: sor_seqsbaij 100 args: -mat_type seqbaij 101 test: 102 suffix: sor_seqsell 103 args: -mat_type seqsell 104 105 TEST*/ 106