1c4762a1bSJed Brown static char help[] = "Demonstrates the use of fast Richardson for SOR. And\n\ 2c4762a1bSJed Brown also tests the MatSOR() routines. Input parameters are:\n\ 3c4762a1bSJed Brown -n <n> : problem dimension\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscksp.h> 6c4762a1bSJed Brown #include <petscpc.h> 7c4762a1bSJed Brown 8d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 9d71ae5a4SJacob Faibussowitsch { 10c4762a1bSJed Brown Mat mat; /* matrix */ 11c4762a1bSJed Brown Vec b, ustar, u; /* vectors (RHS, exact solution, approx solution) */ 12c4762a1bSJed Brown PC pc; /* PC context */ 13c4762a1bSJed Brown KSP ksp; /* KSP context */ 14c4762a1bSJed Brown PetscInt n = 10, i, its, col[3]; 15c4762a1bSJed Brown PetscScalar value[3]; 16c4762a1bSJed Brown KSPType kspname; 17c4762a1bSJed Brown PCType pcname; 18c4762a1bSJed Brown 19327415f7SBarry Smith PetscFunctionBeginUser; 20*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, 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)); 399371c9d4SSatish Balay value[0] = -1.0; 409371c9d4SSatish Balay value[1] = 2.0; 419371c9d4SSatish Balay value[2] = -1.0; 42c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 439371c9d4SSatish Balay col[0] = i - 1; 449371c9d4SSatish Balay col[1] = i; 459371c9d4SSatish Balay col[2] = i + 1; 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 3, col, value, INSERT_VALUES)); 47c4762a1bSJed Brown } 489371c9d4SSatish Balay i = n - 1; 499371c9d4SSatish Balay col[0] = n - 2; 509371c9d4SSatish Balay col[1] = n - 1; 519566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES)); 529371c9d4SSatish Balay i = 0; 539371c9d4SSatish Balay col[0] = 0; 549371c9d4SSatish Balay col[1] = 1; 559371c9d4SSatish Balay value[0] = 2.0; 569371c9d4SSatish Balay value[1] = -1.0; 579566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Compute right-hand-side vector */ 629566063dSJacob Faibussowitsch PetscCall(MatMult(mat, ustar, b)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Create PC context and set up data structures */ 659566063dSJacob Faibussowitsch PetscCall(PCCreate(PETSC_COMM_WORLD, &pc)); 669566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 679566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); 689566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc, mat, mat)); 699566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Create KSP context and set up data structures */ 729566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 739566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPRICHARDSON)); 749566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 759566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc, mat, mat)); 769566063dSJacob Faibussowitsch PetscCall(KSPSetPC(ksp, pc)); 779566063dSJacob Faibussowitsch PetscCall(KSPSetUp(ksp)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Solve the problem */ 809566063dSJacob Faibussowitsch PetscCall(KSPGetType(ksp, &kspname)); 819566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &pcname)); 829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Running %s with %s preconditioning\n", kspname, pcname)); 839566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, b, u)); 849566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &its)); 8563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of iterations %" PetscInt_FMT "\n", its)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Free data structures */ 889566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ustar)); 919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 939566063dSJacob Faibussowitsch PetscCall(PCDestroy(&pc)); 949566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 95b122ec5aSJacob Faibussowitsch return 0; 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown /*TEST 99c4762a1bSJed Brown 100c4762a1bSJed Brown testset: 101c4762a1bSJed Brown args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric 102c4762a1bSJed Brown output_file: output/ex3_1.out 103c4762a1bSJed Brown test: 104c4762a1bSJed Brown suffix: sor_aij 105c4762a1bSJed Brown test: 106c4762a1bSJed Brown suffix: sor_seqbaij 107c4762a1bSJed Brown args: -mat_type seqbaij 108c4762a1bSJed Brown test: 109c4762a1bSJed Brown suffix: sor_seqsbaij 110c4762a1bSJed Brown args: -mat_type seqbaij 111c4762a1bSJed Brown test: 112c4762a1bSJed Brown suffix: sor_seqsell 113c4762a1bSJed Brown args: -mat_type seqsell 114c4762a1bSJed Brown 115c4762a1bSJed Brown TEST*/ 116