Actual source code: ex10.c

  1: /*$Id: ex10.c,v 1.58 2001/09/11 16:33:29 bsmith 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:  #include src/vec/is/impls/general/general.h

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

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

 47:   PetscOptionsHasName(PETSC_NULL,"-table",&table);
 48:   PetscOptionsHasName(PETSC_NULL,"-trans",&trans);
 49:   PetscOptionsHasName(PETSC_NULL,"-partition",&partition);

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

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

 79:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
 80:                            Load system
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 89:     /*
 90:        Load the matrix and vector; then destroy the viewer.
 91:     */
 92:     ierr  = MatLoad(fd,MATMPIAIJ,&A);
 93:     ierr  = PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
 94:     ierrp = VecLoad(fd,&b);
 95:     ierr  = PetscPopErrorHandler();
 96:     if (ierrp) { /* if file contains no RHS, then use a vector of all ones */
 97:       int    m;
 98:       PetscScalar one = 1.0;
 99:       MatGetLocalSize(A,&m,PETSC_NULL);
100:       VecCreate(PETSC_COMM_WORLD,&b);
101:       VecSetSizes(b,m,PETSC_DECIDE);
102:       VecSetFromOptions(b);
103:       VecSet(&one,b);
104:     }
105:     PetscViewerDestroy(fd);

107:     /* 
108:        If the loaded matrix is larger than the vector (due to being padded 
109:        to match the block size of the system), then create a new padded vector.
110:     */
111:     {
112:       int         m,n,j,mvec,start,end,index;
113:       Vec         tmp;
114:       PetscScalar *bold;

116:       /* Create a new vector b by padding the old one */
117:       MatGetLocalSize(A,&m,&n);
118:       VecCreate(PETSC_COMM_WORLD,&tmp);
119:       VecSetSizes(tmp,m,PETSC_DECIDE);
120:       VecSetFromOptions(tmp);
121:       VecGetOwnershipRange(b,&start,&end);
122:       VecGetLocalSize(b,&mvec);
123:       VecGetArray(b,&bold);
124:       for (j=0; j<mvec; j++) {
125:         index = start+j;
126:         ierr  = VecSetValues(tmp,1,&index,bold+j,INSERT_VALUES);
127:       }
128:       VecRestoreArray(b,&bold);
129:       VecDestroy(b);
130:       VecAssemblyBegin(tmp);
131:       VecAssemblyEnd(tmp);
132:       b = tmp;
133:     }
134:     VecDuplicate(b,&x);
135:     VecDuplicate(b,&u);
136:     VecSet(&zero,x);

138:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
139:                       Setup solve for system
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

142:     /*
143:        Conclude profiling last stage; begin profiling next stage.
144:     */
145:     PreLoadStage("SLESSetUp");

147:     if (partition) {
148:       MatPartitioning mpart;
149:       IS              mis,nis,isn,is;
150:       int             *count,size,rank;
151:       Mat             B;
152:       MPI_Comm_size(PETSC_COMM_WORLD,&size);
153:       MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
154:       PetscMalloc(size*sizeof(int),&count);
155:       MatPartitioningCreate(PETSC_COMM_WORLD, &mpart);
156:       MatPartitioningSetAdjacency(mpart, A);
157:       /* MatPartitioningSetVertexWeights(mpart, weight); */
158:       MatPartitioningSetFromOptions(mpart);
159:       MatPartitioningApply(mpart, &mis);
160:       MatPartitioningDestroy(mpart);
161:       ISPartitioningToNumbering(mis,&nis);
162:       ISPartitioningCount(mis,count);
163:       ISDestroy(mis);
164:       ISInvertPermutation(nis, count[rank], &is);
165:       PetscFree(count);
166:       ISDestroy(nis);
167:       ISSort(is);
168:       ISAllGather(is,&isn);
169:       MatGetSubMatrix(A,is,isn,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B);

171:       /* need to move the vector also */
172:       ISDestroy(is);
173:       ISDestroy(isn);
174:       MatDestroy(A);
175:       A    = B;
176:     }
177: 
178:     /*
179:        We also explicitly time this stage via PetscGetTime()
180:     */
181:     PetscGetTime(&tsetup1);

183:     /*
184:        Create linear solver; set operators; set runtime options.
185:     */
186:     SLESCreate(PETSC_COMM_WORLD,&sles);
187:     SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
188:     SLESSetFromOptions(sles);

190:     /* 
191:        Here we explicitly call SLESSetUp() and SLESSetUpOnBlocks() to
192:        enable more precise profiling of setting up the preconditioner.
193:        These calls are optional, since both will be called within
194:        SLESSolve() if they haven't been called already.
195:     */
196:     SLESSetUp(sles,b,x);
197:     SLESSetUpOnBlocks(sles);
198:     PetscGetTime(&tsetup2);
199:     tsetup = tsetup2 - tsetup1;

201:     /*
202:        Tests "diagonal-scaling of preconditioned residual norm" as used 
203:        by many ODE integrator codes including PVODE. Note this is different
204:        than diagonally scaling the matrix before computing the preconditioner
205:     */
206:     PetscOptionsHasName(PETSC_NULL,"-diagonal_scale",&diagonalscale);
207:     if (diagonalscale) {
208:       PC     pc;
209:       int    j,start,end,n;
210:       Vec    scale;
211: 
212:       SLESGetPC(sles,&pc);
213:       VecGetSize(x,&n);
214:       VecDuplicate(x,&scale);
215:       VecGetOwnershipRange(scale,&start,&end);
216:       for (j=start; j<end; j++) {
217:         VecSetValue(scale,j,((PetscReal)(j+1))/((PetscReal)n),INSERT_VALUES);
218:       }
219:       VecAssemblyBegin(scale);
220:       VecAssemblyEnd(scale);
221:       PCDiagonalScaleSet(pc,scale);
222:       VecDestroy(scale);

224:     }

226:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
227:                            Solve system
228:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

230:     /*
231:        Begin profiling next stage
232:     */
233:     PreLoadStage("SLESSolve");

235:     /*
236:        Solve linear system; we also explicitly time this stage.
237:     */
238:     PetscGetTime(&tsolve1);
239:     if (trans) {
240:       SLESSolveTranspose(sles,b,x,&its);
241:     } else {
242:       SLESSolve(sles,b,x,&its);
243:     }
244:     PetscGetTime(&tsolve2);
245:     tsolve = tsolve2 - tsolve1;

247:    /* 
248:        Conclude profiling this stage
249:     */
250:     PreLoadStage("Cleanup");

252:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
253:             Check error, print output, free data structures.
254:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

256:     /* 
257:        Check error
258:     */
259:     if (trans) {
260:       MatMultTranspose(A,x,u);
261:     } else {
262:       MatMult(A,x,u);
263:     }
264:     VecAXPY(&none,b,u);
265:     VecNorm(u,NORM_2,&norm);

267:     /*
268:        Write output (optinally using table for solver details).
269:         - PetscPrintf() handles output for multiprocessor jobs 
270:           by printing from only one processor in the communicator.
271:         - SLESView() prints information about the linear solver.
272:     */
273:     if (table) {
274:       char        *matrixname,slesinfo[120];
275:       PetscViewer viewer;

277:       /*
278:          Open a string viewer; then write info to it.
279:       */
280:       PetscViewerStringOpen(PETSC_COMM_WORLD,slesinfo,120,&viewer);
281:       SLESView(sles,viewer);
282:       PetscStrrchr(file[PreLoadIt],'/',&matrixname);
283:       PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3d %2.0e %2.1e %2.1e %2.1e %s n",
284:                 matrixname,its,norm,tsetup+tsolve,tsetup,tsolve,slesinfo);

286:       /*
287:          Destroy the viewer
288:       */
289:       PetscViewerDestroy(viewer);
290:     } else {
291:       PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
292:       PetscPrintf(PETSC_COMM_WORLD,"Residual norm %An",norm);
293:     }
294:     /* 
295:        Free work space.  All PETSc objects should be destroyed when they
296:        are no longer needed.
297:     */
298:     MatDestroy(A); VecDestroy(b);
299:     VecDestroy(u); VecDestroy(x);
300:     SLESDestroy(sles);
301:   PreLoadEnd();
302:   /* -----------------------------------------------------------
303:                       End of linear solver loop
304:      ----------------------------------------------------------- */

306:   PetscFinalize();
307:   return 0;
308: }