Actual source code: ex10.c

  1: /*$Id: ex10.c,v 1.51 2001/03/23 23:23:55 balay Exp $*/

  3: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.n
  4: This version first preloads and solves a small system, then loads n
  5: another (larger) system and solves it as well.  This example illustratesn
  6: preloading of instructions with the smaller system so that more accuraten
  7: performance monitoring can be done with the larger one (that actuallyn
  8: is the system of interest).  See the 'Performance Hints' chapter of then
  9: users manual for a discussion of preloading.  Input parameters includen
 10:   -f0 <input_file> : first file to load (small system)n
 11:   -f1 <input_file> : second file to load (larger system)nn
 12:   -trans  : solve transpose system insteadnn";

 14: /*T
 15:    Concepts: SLES^solving a linear system
 16:    Processors: n
 17: T*/

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

 29: int main(int argc,char **args)
 30: {
 31:   SLES           sles;             /* linear solver context */
 32:   Mat            A;                /* matrix */
 33:   Vec            x,b,u;          /* approx solution, RHS, exact solution */
 34:   PetscViewer    fd;               /* viewer */
 35:   char           file[2][128];     /* input file name */
 36:   PetscTruth     table,flg,trans;
 37:   int            ierr,its;
 38:   double         norm;
 39:   PetscLogDouble tsetup,tsetup1,tsetup2,tsolve,tsolve1,tsolve2;
 40:   Scalar         zero = 0.0,none = -1.0;
 41:   PetscTruth     preload = PETSC_TRUE,diagonalscale;

 43:   PetscInitialize(&argc,&args,(char *)0,help);

 45:   PetscOptionsHasName(PETSC_NULL,"-table",&table);
 46:   PetscOptionsHasName(PETSC_NULL,"-trans",&trans);

 48:   /* 
 49:      Determine files from which we read the two linear systems
 50:      (matrix and right-hand-side vector).
 51:   */
 52:   PetscOptionsGetString(PETSC_NULL,"-f",file[0],127,&flg);
 53:   if (flg) {
 54:     PetscStrcpy(file[1],file[0]);
 55:   } else {
 56:     PetscOptionsGetString(PETSC_NULL,"-f0",file[0],127,&flg);
 57:     if (!flg) SETERRQ(1,"Must indicate binary file with the -f0 or -f option");
 58:     PetscOptionsGetString(PETSC_NULL,"-f1",file[1],127,&flg);
 59:     if (!flg) {preload = PETSC_FALSE;} /* don't bother with second system */
 60:   }

 62:   /* -----------------------------------------------------------
 63:                   Beginning of linear solver loop
 64:      ----------------------------------------------------------- */
 65:   /* 
 66:      Loop through the linear solve 2 times.  
 67:       - The intention here is to preload and solve a small system;
 68:         then load another (larger) system and solve it as well.
 69:         This process preloads the instructions with the smaller
 70:         system so that more accurate performance monitoring (via
 71:         -log_summary) can be done with the larger one (that actually
 72:         is the system of interest). 
 73:   */
 74:   PreLoadBegin(preload,"Load system");

 76:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
 77:                            Load system
 78:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 80:     /* 
 81:        Open binary file.  Note that we use PETSC_BINARY_RDONLY to indicate
 82:        reading from this file.
 83:     */
 84:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt],PETSC_BINARY_RDONLY,&fd);

 86:     /*
 87:        Load the matrix and vector; then destroy the viewer.
 88:     */
 89:     MatLoad(fd,MATMPIAIJ,&A);
 90:     VecLoad(fd,&b);
 91:     if (ierr) { /* if file contains no RHS, then use a vector of all ones */
 92:       int    m;
 93:       Scalar one = 1.0;
 94:       MatGetLocalSize(A,&m,PETSC_NULL);
 95:       VecCreateMPI(PETSC_COMM_WORLD,m,PETSC_DECIDE,&b);
 96:       VecSet(&one,b);
 97:     }
 98:     PetscViewerDestroy(fd);

100:     /* 
101:        If the loaded matrix is larger than the vector (due to being padded 
102:        to match the block size of the system), then create a new padded vector.
103:     */
104:     {
105:       int    m,n,j,mvec,start,end,index;
106:       Vec    tmp;
107:       Scalar *bold;

109:       /* Create a new vector b by padding the old one */
110:       MatGetLocalSize(A,&m,&n);
111:       VecCreateMPI(PETSC_COMM_WORLD,m,PETSC_DECIDE,&tmp);
112:       VecGetOwnershipRange(b,&start,&end);
113:       VecGetLocalSize(b,&mvec);
114:       VecGetArray(b,&bold);
115:       for (j=0; j<mvec; j++) {
116:         index = start+j;
117:         ierr  = VecSetValues(tmp,1,&index,bold+j,INSERT_VALUES);
118:       }
119:       VecRestoreArray(b,&bold);
120:       VecDestroy(b);
121:       VecAssemblyBegin(tmp);
122:       VecAssemblyEnd(tmp);
123:       b = tmp;
124:     }
125:     VecDuplicate(b,&x);
126:     VecDuplicate(b,&u);
127:     VecSet(&zero,x);

129:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
130:                       Setup solve for system
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

133:     /*
134:        Conclude profiling last stage; begin profiling next stage.
135:     */
136:     PreLoadStage("SLESSetUp");

138:     /*
139:        We also explicitly time this stage via PetscGetTime()
140:     */
141:     PetscGetTime(&tsetup1);

143:     /*
144:        Create linear solver; set operators; set runtime options.
145:     */
146:     SLESCreate(PETSC_COMM_WORLD,&sles);
147:     SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
148:     SLESSetFromOptions(sles);

150:     /* 
151:        Here we explicitly call SLESSetUp() and SLESSetUpOnBlocks() to
152:        enable more precise profiling of setting up the preconditioner.
153:        These calls are optional, since both will be called within
154:        SLESSolve() if they haven't been called already.
155:     */
156:     SLESSetUp(sles,b,x);
157:     SLESSetUpOnBlocks(sles);
158:     PetscGetTime(&tsetup2);
159:     tsetup = tsetup2 - tsetup1;

161:     /*
162:        Tests "diagonal-scaling of preconditioned residual norm" as used 
163:        by many ODE integrator codes including PVODE. Note this is different
164:        than diagonally scaling the matrix before computing the preconditioner
165:     */
166:     PetscOptionsHasName(PETSC_NULL,"-diagonal_scale",&diagonalscale);
167:     if (diagonalscale) {
168:       PC     pc;
169:       int    j,start,end,n;
170:       Vec    scale;
171: 
172:       SLESGetPC(sles,&pc);
173:       VecGetSize(x,&n);
174:       VecDuplicate(x,&scale);
175:       VecGetOwnershipRange(scale,&start,&end);
176:       for (j=start; j<end; j++) {
177:         VecSetValue(scale,j,((double)(j+1))/((double)n),INSERT_VALUES);
178:       }
179:       VecAssemblyBegin(scale);
180:       VecAssemblyEnd(scale);
181:       PCDiagonalScaleSet(pc,scale);
182:       VecDestroy(scale);

184:     }

186:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
187:                            Solve system
188:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

190:     /*
191:        Begin profiling next stage
192:     */
193:     PreLoadStage("SLESSolve");

195:     /*
196:        Solve linear system; we also explicitly time this stage.
197:     */
198:     PetscGetTime(&tsolve1);
199:     if (trans) {
200:       SLESSolveTranspose(sles,b,x,&its);
201:     } else {
202:       SLESSolve(sles,b,x,&its);
203:     }
204:     PetscGetTime(&tsolve2);
205:     tsolve = tsolve2 - tsolve1;

207:    /* 
208:        Conclude profiling this stage
209:     */
210:     PreLoadStage("Cleanup");

212:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
213:             Check error, print output, free data structures.
214:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

216:     /* 
217:        Check error
218:     */
219:     if (trans) {
220:       MatMultTranspose(A,x,u);
221:     } else {
222:       MatMult(A,x,u);
223:     }
224:     VecAXPY(&none,b,u);
225:     VecNorm(u,NORM_2,&norm);

227:     /*
228:        Write output (optinally using table for solver details).
229:         - PetscPrintf() handles output for multiprocessor jobs 
230:           by printing from only one processor in the communicator.
231:         - SLESView() prints information about the linear solver.
232:     */
233:     if (table) {
234:       char        *matrixname,slesinfo[120];
235:       PetscViewer viewer;

237:       /*
238:          Open a string viewer; then write info to it.
239:       */
240:       PetscViewerStringOpen(PETSC_COMM_WORLD,slesinfo,120,&viewer);
241:       SLESView(sles,viewer);
242:       PetscStrrchr(file[PreLoadIt],'/',&matrixname);
243:       PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3d %2.0e %2.1e %2.1e %2.1e %s n",
244:                 matrixname,its,norm,tsetup+tsolve,tsetup,tsolve,slesinfo);

246:       /*
247:          Destroy the viewer
248:       */
249:       PetscViewerDestroy(viewer);
250:     } else {
251:       PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
252:       PetscPrintf(PETSC_COMM_WORLD,"Residual norm %An",norm);
253:     }
254:     /* 
255:        Free work space.  All PETSc objects should be destroyed when they
256:        are no longer needed.
257:     */
258:     MatDestroy(A); VecDestroy(b);
259:     VecDestroy(u); VecDestroy(x);
260:     SLESDestroy(sles);
261:   PreLoadEnd();
262:   /* -----------------------------------------------------------
263:                       End of linear solver loop
264:      ----------------------------------------------------------- */

266:   PetscFinalize();
267:   return 0;
268: }