Actual source code: ex48.c


  2: static char help[] = "Solves a tridiagonal linear system with KSP.\n\n";

  4: /*
  5:     Test example that demonstrates how MINRES can produce a dp of zero
  6:     but still converge.

  8:     Provided by: Mark Filipiak <mjf@staffmail.ed.ac.uk>
  9: */
 10: #include <petscksp.h>

 12: int main(int argc,char **args)
 13: {
 14:   Vec            x, b, u;      /* approx solution, RHS, exact solution */
 15:   Mat            A;            /* linear system matrix */
 16:   KSP            ksp;         /* linear solver context */
 17:   PC             pc;           /* preconditioner context */
 18:   PetscReal      norm;
 19:   PetscInt       i,n = 2,col[3],its;
 20:   PetscMPIInt    size;
 21:   PetscScalar    one = 1.0,value[3];
 22:   PetscBool      nonzeroguess = PETSC_FALSE;

 24:   PetscInitialize(&argc,&args,(char*)0,help);
 25:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   PetscOptionsGetBool(NULL,NULL,"-nonzero_guess",&nonzeroguess,NULL);

 30:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 31:          Compute the matrix and right-hand-side vector that define
 32:          the linear system, Ax = b.
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 35:   /*
 36:      Create vectors.  Note that we form 1 vector from scratch and
 37:      then duplicate as needed.
 38:   */
 39:   VecCreate(PETSC_COMM_WORLD,&x);
 40:   PetscObjectSetName((PetscObject) x, "Solution");
 41:   VecSetSizes(x,PETSC_DECIDE,n);
 42:   VecSetFromOptions(x);
 43:   VecDuplicate(x,&b);
 44:   VecDuplicate(x,&u);

 46:   /*
 47:      Create matrix.  When using MatCreate(), the matrix format can
 48:      be specified at runtime.

 50:      Performance tuning note:  For problems of substantial size,
 51:      preallocation of matrix memory is crucial for attaining good
 52:      performance. See the matrix chapter of the users manual for details.
 53:   */
 54:   MatCreate(PETSC_COMM_WORLD,&A);
 55:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 56:   MatSetFromOptions(A);
 57:   MatSetUp(A);

 59:   /*
 60:      Assemble matrix
 61:   */
 62:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 63:   for (i=1; i<n-1; i++) {
 64:     col[0] = i-1; col[1] = i; col[2] = i+1;
 65:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 66:   }
 67:   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
 68:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 69:   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 70:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 71:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 72:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 74:   /*
 75:      Set constant right-hand-side vector.
 76:   */
 77:   VecSet(b,one);
 78:   /*
 79:      Solution = RHS for the matrix [[2 -1] [-1 2]] and RHS [1 1]
 80:   */
 81:   VecSet(u,one);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                 Create the linear solver and set various options
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 86:   /*
 87:      Create linear solver context
 88:   */
 89:   KSPCreate(PETSC_COMM_WORLD,&ksp);

 91:   /*
 92:      Set operators. Here the matrix that defines the linear system
 93:      also serves as the preconditioning matrix.
 94:   */
 95:   KSPSetOperators(ksp,A,A);

 97:   /*
 98:      Set linear solver defaults for this problem (optional).
 99:      - By extracting the KSP and PC contexts from the KSP context,
100:        we can then directly call any KSP and PC routines to set
101:        various options.
102:      - The following four statements are optional; all of these
103:        parameters could alternatively be specified at runtime via
104:        KSPSetFromOptions();
105:   */
106:   KSPGetPC(ksp,&pc);
107:   PCSetType(pc,PCJACOBI);
108:   KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

110:   /*
111:     Set runtime options, e.g.,
112:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
113:     These options will override those specified above as long as
114:     KSPSetFromOptions() is called _after_ any other customization
115:     routines.
116:   */
117:   KSPSetFromOptions(ksp);

119:   if (nonzeroguess) {
120:     PetscScalar p = .5;
121:     VecSet(x,p);
122:     KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
123:   }

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126:                       Solve the linear system
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128:   /*
129:      Solve linear system
130:   */
131:   KSPSolve(ksp,b,x);

133:   /*
134:      View solver info; we could instead use the option -ksp_view to
135:      print this info to the screen at the conclusion of KSPSolve().
136:   */
137:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:                       Check solution and clean up
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142:   /*
143:      Check the error
144:   */
145:   VecAXPY(x,-1.0,u);
146:   VecNorm(x,NORM_2,&norm);
147:   KSPGetIterationNumber(ksp,&its);
148:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);

150:   /*
151:      Free work space.  All PETSc objects should be destroyed when they
152:      are no longer needed.
153:   */
154:   VecDestroy(&x)); PetscCall(VecDestroy(&u);
155:   VecDestroy(&b)); PetscCall(MatDestroy(&A);
156:   KSPDestroy(&ksp);

158:   /*
159:      Always call PetscFinalize() before exiting a program.  This routine
160:        - finalizes the PETSc libraries as well as MPI
161:        - provides summary and diagnostic information if certain runtime
162:          options are chosen (e.g., -log_view).
163:   */
164:   PetscFinalize();
165:   return 0;
166: }

168: /*TEST

170:    test:
171:       args: -ksp_monitor_short -ksp_converged_reason -ksp_type minres -pc_type jacobi -ksp_error_if_not_converged

173: TEST*/