Actual source code: ex48.c

petsc-3.7.4 2016-10-02
Report Typos and Errors
  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>
 11: #include <petsc/private/kspimpl.h>

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

 28:   PetscInitialize(&argc,&args,(char*)0,help);
 29:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 30:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only!");
 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetBool(NULL,NULL,"-nonzero_guess",&nonzeroguess,NULL);


 35:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 36:          Compute the matrix and right-hand-side vector that define
 37:          the linear system, Ax = b.
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 51:   /*
 52:      Create matrix.  When using MatCreate(), the matrix format can
 53:      be specified at runtime.

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

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

 79:   /*
 80:      Set constant right-hand-side vector.
 81:   */
 82:   VecSet(b,one);
 83:   /*
 84:      Solution = RHS for the matrix [[2 -1] [-1 2]] and RHS [1 1]
 85:   */
 86:   VecSet(u,one);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                 Create the linear solver and set various options
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 91:   /*
 92:      Create linear solver context
 93:   */
 94:   KSPCreate(PETSC_COMM_WORLD,&ksp);

 96:   /*
 97:      Set operators. Here the matrix that defines the linear system
 98:      also serves as the preconditioning matrix.
 99:   */
100:   KSPSetOperators(ksp,A,A);

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

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

124:   if (nonzeroguess) {
125:     PetscScalar p = .5;
126:     VecSet(x,p);
127:     KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
128:   }

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:                       Solve the linear system
132:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133:   /*
134:      Solve linear system
135:   */
136:   KSPSolve(ksp,b,x);

138:   /*
139:      View solver info; we could instead use the option -ksp_view to
140:      print this info to the screen at the conclusion of KSPSolve().
141:   */
142:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:                       Check solution and clean up
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147:   /*
148:      Check the error
149:   */
150:   VecAXPY(x,neg_one,u);
151:   VecNorm(x,NORM_2,&norm);
152:   KSPGetIterationNumber(ksp,&its);
153:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);


156:   /*
157:      Free work space.  All PETSc objects should be destroyed when they
158:      are no longer needed.
159:   */
160:   VecDestroy(&x); VecDestroy(&u);
161:   VecDestroy(&b); MatDestroy(&A);
162:   KSPDestroy(&ksp);

164:   /*
165:      Always call PetscFinalize() before exiting a program.  This routine
166:        - finalizes the PETSc libraries as well as MPI
167:        - provides summary and diagnostic information if certain runtime
168:          options are chosen (e.g., -log_summary).
169:   */
170:   PetscFinalize();
171:   return 0;
172: }