Actual source code: ex5.c

  1: /*$Id: ex5.c,v 1.93 2001/09/11 16:33:29 bsmith Exp $*/

  3: static char help[] = "Solves two linear systems in parallel with SLES.  The coden
  4: illustrates repeated solution of linear systems with the same preconditionern
  5: method but different matrices (having the same nonzero structure).  The coden
  6: also uses multiple profiling stages.  Input arguments aren
  7:   -m <size> : problem sizen
  8:   -mat_nonsym : use nonsymmetric matrix (default is symmetric)nn";

 10: /*T
 11:    Concepts: SLES^repeatedly solving linear systems;
 12:    Concepts: PetscLog^profiling multiple stages of code;
 13:    Processors: n
 14: T*/

 16: /* 
 17:   Include "petscsles.h" so that we can use SLES solvers.  Note that this file
 18:   automatically includes:
 19:      petsc.h       - base PETSc routines   petscvec.h - vectors
 20:      petscsys.h    - system routines       petscmat.h - matrices
 21:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 22:      petscviewer.h - viewers               petscpc.h  - preconditioners
 23: */
 24:  #include petscsles.h

 26: int main(int argc,char **args)
 27: {
 28:   SLES         sles;             /* linear solver context */
 29:   Mat          C;                /* matrix */
 30:   Vec          x,u,b;          /* approx solution, RHS, exact solution */
 31:   PetscReal    norm;             /* norm of solution error */
 32:   PetscScalar  v,none = -1.0;
 33:   int          I,J,ldim,ierr,low,high,iglobal,Istart,Iend;
 34:   int          i,j,m = 3,n = 2,rank,size,its;
 35:   PetscTruth   mat_nonsymmetric;
 36:   int          stages[2];

 38:   PetscInitialize(&argc,&args,(char *)0,help);
 39:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 40:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 41:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 42:   n = 2*size;

 44:   /*
 45:      Set flag if we are doing a nonsymmetric problem; the default is symmetric.
 46:   */
 47:   PetscOptionsHasName(PETSC_NULL,"-mat_nonsym",&mat_nonsymmetric);

 49:   /*
 50:      Register two stages for separate profiling of the two linear solves.
 51:      Use the runtime option -log_summary for a printout of performance
 52:      statistics at the program's conlusion.
 53:   */
 54:   PetscLogStageRegister(&stages[0],"Original Solve");
 55:   PetscLogStageRegister(&stages[1],"Second Solve");

 57:   /* -------------- Stage 0: Solve Original System ---------------------- */
 58:   /* 
 59:      Indicate to PETSc profiling that we're beginning the first stage
 60:   */
 61:   PetscLogStagePush(stages[0]);

 63:   /* 
 64:      Create parallel matrix, specifying only its global dimensions.
 65:      When using MatCreate(), the matrix format can be specified at
 66:      runtime. Also, the parallel partitioning of the matrix is
 67:      determined by PETSc at runtime.
 68:   */
 69:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
 70:   MatSetFromOptions(C);

 72:   /* 
 73:      Currently, all PETSc parallel matrix formats are partitioned by
 74:      contiguous chunks of rows across the processors.  Determine which
 75:      rows of the matrix are locally owned. 
 76:   */
 77:   MatGetOwnershipRange(C,&Istart,&Iend);

 79:   /* 
 80:      Set matrix entries matrix in parallel.
 81:       - Each processor needs to insert only elements that it owns
 82:         locally (but any non-local elements will be sent to the
 83:         appropriate processor during matrix assembly). 
 84:       - Always specify global row and columns of matrix entries.
 85:   */
 86:   for (I=Istart; I<Iend; I++) {
 87:     v = -1.0; i = I/n; j = I - i*n;
 88:     if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 89:     if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 90:     if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 91:     if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 92:     v = 4.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
 93:   }

 95:   /*
 96:      Make the matrix nonsymmetric if desired
 97:   */
 98:   if (mat_nonsymmetric) {
 99:     for (I=Istart; I<Iend; I++) {
100:       v = -1.5; i = I/n;
101:       if (i>1)   {J = I-n-1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
102:     }
103:   } else {
104:     MatSetOption(C,MAT_SYMMETRIC);
105:   }

107:   /* 
108:      Assemble matrix, using the 2-step process:
109:        MatAssemblyBegin(), MatAssemblyEnd()
110:      Computations can be done while messages are in transition
111:      by placing code between these two statements.
112:   */
113:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
114:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

116:   /* 
117:      Create parallel vectors.
118:       - When using VecSetSizes(), we specify only the vector's global
119:         dimension; the parallel partitioning is determined at runtime. 
120:       - Note: We form 1 vector from scratch and then duplicate as needed.
121:   */
122:   VecCreate(PETSC_COMM_WORLD,&u);
123:   VecSetSizes(u,PETSC_DECIDE,m*n);
124:   VecSetFromOptions(u);
125:   VecDuplicate(u,&b);
126:   VecDuplicate(b,&x);

128:   /* 
129:      Currently, all parallel PETSc vectors are partitioned by
130:      contiguous chunks across the processors.  Determine which
131:      range of entries are locally owned.
132:   */
133:   VecGetOwnershipRange(x,&low,&high);

135:   /*
136:     Set elements within the exact solution vector in parallel.
137:      - Each processor needs to insert only elements that it owns
138:        locally (but any non-local entries will be sent to the
139:        appropriate processor during vector assembly).
140:      - Always specify global locations of vector entries.
141:   */
142:   VecGetLocalSize(x,&ldim);
143:   for (i=0; i<ldim; i++) {
144:     iglobal = i + low;
145:     v = (PetscScalar)(i + 100*rank);
146:     VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);
147:   }

149:   /* 
150:      Assemble vector, using the 2-step process:
151:        VecAssemblyBegin(), VecAssemblyEnd()
152:      Computations can be done while messages are in transition,
153:      by placing code between these two statements.
154:   */
155:   VecAssemblyBegin(u);
156:   VecAssemblyEnd(u);

158:   /* 
159:      Compute right-hand-side vector
160:   */
161:   MatMult(C,u,b);
162: 
163:   /* 
164:     Create linear solver context
165:   */
166:   SLESCreate(PETSC_COMM_WORLD,&sles);

168:   /* 
169:      Set operators. Here the matrix that defines the linear system
170:      also serves as the preconditioning matrix.
171:   */
172:   SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);

174:   /* 
175:      Set runtime options (e.g., -ksp_type <type> -pc_type <type>)
176:   */

178:   SLESSetFromOptions(sles);

180:   /* 
181:      Solve linear system.  Here we explicitly call SLESSetUp() for more
182:      detailed performance monitoring of certain preconditioners, such
183:      as ICC and ILU.  This call is optional, as SLESSetUp() will
184:      automatically be called within SLESSolve() if it hasn't been
185:      called already.
186:   */
187:   SLESSetUp(sles,b,x);
188:   SLESSolve(sles,b,x,&its);
189: 
190:   /* 
191:      Check the error
192:   */
193:   VecAXPY(&none,u,x);
194:   VecNorm(x,NORM_2,&norm);
195:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %dn",norm,its);

197:   /* -------------- Stage 1: Solve Second System ---------------------- */
198:   /* 
199:      Solve another linear system with the same method.  We reuse the SLES
200:      context, matrix and vector data structures, and hence save the
201:      overhead of creating new ones.

203:      Indicate to PETSc profiling that we're concluding the first
204:      stage with PetscLogStagePop(), and beginning the second stage with
205:      PetscLogStagePush().
206:   */
207:   PetscLogStagePop();
208:   PetscLogStagePush(stages[1]);

210:   /* 
211:      Initialize all matrix entries to zero.  MatZeroEntries() retains the
212:      nonzero structure of the matrix for sparse formats.
213:   */
214:   MatZeroEntries(C);

216:   /* 
217:      Assemble matrix again.  Note that we retain the same matrix data
218:      structure and the same nonzero pattern; we just change the values
219:      of the matrix entries.
220:   */
221:   for (i=0; i<m; i++) {
222:     for (j=2*rank; j<2*rank+2; j++) {
223:       v = -1.0;  I = j + n*i;
224:       if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
225:       if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
226:       if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
227:       if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
228:       v = 6.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
229:     }
230:   }
231:   if (mat_nonsymmetric) {
232:     for (I=Istart; I<Iend; I++) {
233:       v = -1.5; i = I/n;
234:       if (i>1)   {J = I-n-1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
235:     }
236:   }
237:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
238:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

240:   /* 
241:      Compute another right-hand-side vector
242:   */
243:   MatMult(C,u,b);

245:   /* 
246:      Set operators. Here the matrix that defines the linear system
247:      also serves as the preconditioning matrix.
248:       - The flag SAME_NONZERO_PATTERN indicates that the
249:         preconditioning matrix has identical nonzero structure
250:         as during the last linear solve (although the values of
251:         the entries have changed). Thus, we can save some
252:         work in setting up the preconditioner (e.g., no need to
253:         redo symbolic factorization for ILU/ICC preconditioners).
254:       - If the nonzero structure of the matrix is different during
255:         the second linear solve, then the flag DIFFERENT_NONZERO_PATTERN
256:         must be used instead.  If you are unsure whether the
257:         matrix structure has changed or not, use the flag
258:         DIFFERENT_NONZERO_PATTERN.
259:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
260:         believes your assertion and does not check the structure
261:         of the matrix.  If you erroneously claim that the structure
262:         is the same when it actually is not, the new preconditioner
263:         will not function correctly.  Thus, use this optimization
264:         feature with caution!
265:   */
266:   SLESSetOperators(sles,C,C,SAME_NONZERO_PATTERN);

268:   /* 
269:      Solve linear system
270:   */
271:   SLESSetUp(sles,b,x);
272:   SLESSolve(sles,b,x,&its);

274:   /* 
275:      Check the error
276:   */
277:   VecAXPY(&none,u,x);
278:   VecNorm(x,NORM_2,&norm);
279:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %dn",norm,its);

281:   /* 
282:      Free work space.  All PETSc objects should be destroyed when they
283:      are no longer needed.
284:   */
285:   SLESDestroy(sles);
286:   VecDestroy(u);
287:   VecDestroy(x);
288:   VecDestroy(b);
289:   MatDestroy(C);

291:   /*
292:      Indicate to PETSc profiling that we're concluding the second stage 
293:   */
294:   PetscLogStagePop();

296:   PetscFinalize();
297:   return 0;
298: }