Actual source code: ex3.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  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";

  6:  #include <petscksp.h>
  7:  #include <petscpc.h>

  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 */
 16:   PetscInt       n = 10,i,its,col[3];
 17:   PetscScalar    value[3];
 18:   KSPType        kspname;
 19:   PCType         pcname;

 21:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 22:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 24:   /* Create and initialize vectors */
 25:   VecCreateSeq(PETSC_COMM_SELF,n,&b);
 26:   VecCreateSeq(PETSC_COMM_SELF,n,&ustar);
 27:   VecCreateSeq(PETSC_COMM_SELF,n,&u);
 28:   VecSet(ustar,1.0);
 29:   VecSet(u,0.0);

 31:   /* Create and assemble matrix */
 32:   MatCreate(PETSC_COMM_SELF,&mat);
 33:   MatSetType(mat,MATSEQAIJ);
 34:   MatSetSizes(mat,n,n,n,n);
 35:   MatSetFromOptions(mat);
 36:   MatSeqAIJSetPreallocation(mat,3,NULL);
 37:   MatSeqBAIJSetPreallocation(mat,1,3,NULL);
 38:   MatSeqSBAIJSetPreallocation(mat,1,3,NULL);
 39:   MatSeqSELLSetPreallocation(mat,3,NULL);
 40:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 41:   for (i=1; i<n-1; i++) {
 42:     col[0] = i-1; col[1] = i; col[2] = i+1;
 43:     MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES);
 44:   }
 45:   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
 46:   MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
 47:   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 48:   MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
 49:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 50:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

 52:   /* Compute right-hand-side vector */
 53:   MatMult(mat,ustar,b);

 55:   /* Create PC context and set up data structures */
 56:   PCCreate(PETSC_COMM_WORLD,&pc);
 57:   PCSetType(pc,PCNONE);
 58:   PCSetFromOptions(pc);
 59:   PCSetOperators(pc,mat,mat);
 60:   PCSetUp(pc);

 62:   /* Create KSP context and set up data structures */
 63:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 64:   KSPSetType(ksp,KSPRICHARDSON);
 65:   KSPSetFromOptions(ksp);
 66:   PCSetOperators(pc,mat,mat);
 67:   KSPSetPC(ksp,pc);
 68:   KSPSetUp(ksp);

 70:   /* Solve the problem */
 71:   KSPGetType(ksp,&kspname);
 72:   PCGetType(pc,&pcname);
 73:   PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname);
 74:   KSPSolve(ksp,b,u);
 75:   KSPGetIterationNumber(ksp,&its);
 76:   PetscPrintf(PETSC_COMM_WORLD,"Number of iterations %D\n",its);

 78:   /* Free data structures */
 79:   KSPDestroy(&ksp);
 80:   VecDestroy(&u);
 81:   VecDestroy(&ustar);
 82:   VecDestroy(&b);
 83:   MatDestroy(&mat);
 84:   PCDestroy(&pc);
 85:   PetscFinalize();
 86:   return ierr;
 87: }

 89: /*TEST

 91:    testset:
 92:      args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric
 93:      output_file: output/ex3_1.out
 94:      test:
 95:        suffix: sor_aij
 96:      test:
 97:        suffix: sor_seqbaij
 98:        args: -mat_type seqbaij
 99:      test:
100:        suffix: sor_seqsbaij
101:        args: -mat_type seqbaij
102:      test:
103:        suffix: sor_seqsell
104:        args: -mat_type seqsell

106: TEST*/