Actual source code: ex10.c

  2: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.\n\
  3: This version first preloads and solves a small system, then loads \n\
  4: another (larger) system and solves it as well.  This example illustrates\n\
  5: preloading of instructions with the smaller system so that more accurate\n\
  6: performance monitoring can be done with the larger one (that actually\n\
  7: is the system of interest).  See the 'Performance Hints' chapter of the\n\
  8: users manual for a discussion of preloading.  Input parameters include\n\
  9:   -f0 <input_file> : first file to load (small system)\n\
 10:   -f1 <input_file> : second file to load (larger system)\n\n\
 11:   -trans  : solve transpose system instead\n\n";
 12: /*
 13:   This code can be used to test PETSc interface to other packages.\n\
 14:   Examples of command line options:       \n\
 15:    ex10 -f0 <datafile> -ksp_type preonly  \n\
 16:         -help -ksp_view                  \n\
 17:         -num_numfac <num_numfac> -num_rhs <num_rhs> \n\
 18:         -ksp_type preonly -pc_type lu -mat_type aijspooles/superlu/superlu_dist/aijmumps \n\
 19:         -ksp_type preonly -pc_type cholesky -mat_type sbaijspooles/dscpack/sbaijmumps \n\
 20:         -f0 <A> -fB <B> -mat_type sbaijmumps -ksp_type preonly -pc_type cholesky -test_inertia -mat_sigma <sigma> \n\
 21:    mpirun -np <np> ex10 -f0 <datafile> -ksp_type cg -pc_type asm -pc_asm_type basic -sub_pc_type icc -mat_type sbaij
 22:  \n\n";
 23: */
 24: /*T
 25:    Concepts: KSP^solving a linear system
 26:    Processors: n
 27: T*/

 29: /* 
 30:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 31:   automatically includes:
 32:      petsc.h       - base PETSc routines   petscvec.h - vectors
 33:      petscsys.h    - system routines       petscmat.h - matrices
 34:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 35:      petscviewer.h - viewers               petscpc.h  - preconditioners
 36: */
 37:  #include petscksp.h


 42: int main(int argc,char **args)
 43: {
 44:   KSP            ksp;             /* linear solver context */
 45:   Mat            A,B;            /* matrix */
 46:   Vec            x,b,u;          /* approx solution, RHS, exact solution */
 47:   PetscViewer    fd;               /* viewer */
 48:   char           file[3][PETSC_MAX_PATH_LEN];     /* input file name */
 49:   PetscTruth     table,flg,flgB=PETSC_FALSE,trans=PETSC_FALSE,partition=PETSC_FALSE;
 50:   PetscErrorCode ierr,ierrp;
 51:   PetscInt       its;
 52:   PetscReal      norm;
 53:   PetscLogDouble tsetup,tsetup1,tsetup2,tsolve,tsolve1,tsolve2;
 54:   PetscScalar    zero = 0.0,none = -1.0;
 55:   PetscTruth     preload=PETSC_TRUE,diagonalscale,isSymmetric,cknorm=PETSC_FALSE,Test_MatDuplicate=PETSC_FALSE;
 56:   PetscInt       num_numfac;
 57:   PetscScalar    sigma;

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

 61:   PetscOptionsHasName(PETSC_NULL,"-table",&table);
 62:   PetscOptionsHasName(PETSC_NULL,"-trans",&trans);
 63:   PetscOptionsHasName(PETSC_NULL,"-partition",&partition);

 65:   /* 
 66:      Determine files from which we read the two linear systems
 67:      (matrix and right-hand-side vector).
 68:   */
 69:   PetscOptionsGetString(PETSC_NULL,"-f",file[0],PETSC_MAX_PATH_LEN-1,&flg);
 70:   if (flg) {
 71:     PetscStrcpy(file[1],file[0]);
 72:     preload = PETSC_FALSE;
 73:   } else {
 74:     PetscOptionsGetString(PETSC_NULL,"-f0",file[0],PETSC_MAX_PATH_LEN-1,&flg);
 75:     if (!flg) SETERRQ(1,"Must indicate binary file with the -f0 or -f option");
 76:     PetscOptionsGetString(PETSC_NULL,"-f1",file[1],PETSC_MAX_PATH_LEN-1,&flg);
 77:     if (!flg) {preload = PETSC_FALSE;} /* don't bother with second system */
 78:   }

 80:   /* -----------------------------------------------------------
 81:                   Beginning of linear solver loop
 82:      ----------------------------------------------------------- */
 83:   /* 
 84:      Loop through the linear solve 2 times.  
 85:       - The intention here is to preload and solve a small system;
 86:         then load another (larger) system and solve it as well.
 87:         This process preloads the instructions with the smaller
 88:         system so that more accurate performance monitoring (via
 89:         -log_summary) can be done with the larger one (that actually
 90:         is the system of interest). 
 91:   */
 92:   PreLoadBegin(preload,"Load system");

 94:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
 95:                            Load system
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:     /* 
 99:        Open binary file.  Note that we use PETSC_FILE_RDONLY to indicate
100:        reading from this file.
101:     */
102:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt],PETSC_FILE_RDONLY,&fd);

104:     /*
105:        Load the matrix and vector; then destroy the viewer.
106:     */
107:     MatLoad(fd,MATAIJ,&A);
108:     PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
109:     ierrp = VecLoad(fd,PETSC_NULL,&b);
110:     PetscPopErrorHandler();
111:     if (ierrp) { /* if file contains no RHS, then use a vector of all ones */
112:       PetscInt    m;
113:       PetscScalar one = 1.0;
114:       PetscLogInfo(0,"Using vector of ones for RHS\n");
115:       MatGetLocalSize(A,&m,PETSC_NULL);
116:       VecCreate(PETSC_COMM_WORLD,&b);
117:       VecSetSizes(b,m,PETSC_DECIDE);
118:       VecSetFromOptions(b);
119:       VecSet(&one,b);
120:     }
121:     PetscViewerDestroy(fd);

123:     /* Test MatDuplicate() */
124:     if (Test_MatDuplicate){
125:       MatDuplicate(A,MAT_COPY_VALUES,&B);
126:       MatEqual(A,B,&flg);
127:       if (!flg){
128:         PetscPrintf(PETSC_COMM_WORLD,"  A != B \n");
129:       }
130:       MatDestroy(B);
131:     }

133:     /* Add a shift to A */
134:     PetscOptionsGetScalar(PETSC_NULL,"-mat_sigma",&sigma,&flg);
135:     if(flg) {
136:       PetscOptionsGetString(PETSC_NULL,"-fB",file[2],PETSC_MAX_PATH_LEN-1,&flgB);
137:       if (flgB){
138:         /* load B to get A = A + sigma*B */
139:         PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],PETSC_FILE_RDONLY,&fd);
140:         MatLoad(fd,MATAIJ,&B);
141:         PetscViewerDestroy(fd);
142:         MatAXPY(&sigma,B,A,DIFFERENT_NONZERO_PATTERN); /* A <- sigma*B + A */
143:       } else {
144:         MatShift(&sigma,A);
145:       }
146:     }

148:     /* Check whether A is symmetric */
149:     PetscOptionsHasName(PETSC_NULL, "-check_symmetry", &flg);
150:     if (flg) {
151:       Mat Atrans;
152:       MatTranspose(A, &Atrans);
153:       MatEqual(A, Atrans, &isSymmetric);
154:       if (isSymmetric) {
155:         PetscPrintf(PETSC_COMM_WORLD,"A is symmetric \n");
156:       } else {
157:         PetscPrintf(PETSC_COMM_WORLD,"A is non-symmetric \n");
158:       }
159:       MatDestroy(Atrans);
160:     }

162:     /* 
163:        If the loaded matrix is larger than the vector (due to being padded 
164:        to match the block size of the system), then create a new padded vector.
165:     */
166:     {
167:       PetscInt    m,n,j,mvec,start,end,indx;
168:       Vec         tmp;
169:       PetscScalar *bold;

171:       /* Create a new vector b by padding the old one */
172:       MatGetLocalSize(A,&m,&n);
173:       VecCreate(PETSC_COMM_WORLD,&tmp);
174:       VecSetSizes(tmp,m,PETSC_DECIDE);
175:       VecSetFromOptions(tmp);
176:       VecGetOwnershipRange(b,&start,&end);
177:       VecGetLocalSize(b,&mvec);
178:       VecGetArray(b,&bold);
179:       for (j=0; j<mvec; j++) {
180:         indx = start+j;
181:         VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
182:       }
183:       VecRestoreArray(b,&bold);
184:       VecDestroy(b);
185:       VecAssemblyBegin(tmp);
186:       VecAssemblyEnd(tmp);
187:       b = tmp;
188:     }
189:     VecDuplicate(b,&x);
190:     VecDuplicate(b,&u);
191:     VecSet(&zero,x);

193:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
194:                       Setup solve for system
195:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */


198:     if (partition) {
199:       MatPartitioning mpart;
200:       IS              mis,nis,isn,is;
201:       PetscInt        *count;
202:       PetscMPIInt     rank,size;
203:       Mat             BB;
204:       MPI_Comm_size(PETSC_COMM_WORLD,&size);
205:       MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
206:       PetscMalloc(size*sizeof(PetscInt),&count);
207:       MatPartitioningCreate(PETSC_COMM_WORLD, &mpart);
208:       MatPartitioningSetAdjacency(mpart, A);
209:       /* MatPartitioningSetVertexWeights(mpart, weight); */
210:       MatPartitioningSetFromOptions(mpart);
211:       MatPartitioningApply(mpart, &mis);
212:       MatPartitioningDestroy(mpart);
213:       ISPartitioningToNumbering(mis,&nis);
214:       ISPartitioningCount(mis,count);
215:       ISDestroy(mis);
216:       ISInvertPermutation(nis, count[rank], &is);
217:       PetscFree(count);
218:       ISDestroy(nis);
219:       ISSort(is);
220:       ISAllGather(is,&isn);
221:       MatGetSubMatrix(A,is,isn,PETSC_DECIDE,MAT_INITIAL_MATRIX,&BB);

223:       /* need to move the vector also */
224:       ISDestroy(is);
225:       ISDestroy(isn);
226:       MatDestroy(A);
227:       A    = BB;
228:     }
229: 
230:     /*
231:        Conclude profiling last stage; begin profiling next stage.
232:     */
233:     PreLoadStage("KSPSetUp");

235:     /*
236:        We also explicitly time this stage via PetscGetTime()
237:     */
238:     PetscGetTime(&tsetup1);

240:     /*
241:        Create linear solver; set operators; set runtime options.
242:     */
243:     KSPCreate(PETSC_COMM_WORLD,&ksp);

245:     num_numfac = 1;
246:     PetscOptionsGetInt(PETSC_NULL,"-num_numfac",&num_numfac,PETSC_NULL);
247:     while ( num_numfac-- ){
248:       /* KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN); */
249:     KSPSetOperators(ksp,A,A,SAME_NONZERO_PATTERN);
250:     KSPSetFromOptions(ksp);

252:     /* 
253:        Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
254:        enable more precise profiling of setting up the preconditioner.
255:        These calls are optional, since both will be called within
256:        KSPSolve() if they haven't been called already.
257:     */
258:     KSPSetUp(ksp);
259:     KSPSetUpOnBlocks(ksp);
260:     PetscGetTime(&tsetup2);
261:     tsetup = tsetup2 - tsetup1;

263:     /*
264:       Test MatGetInertia()
265:       Usage:
266:       ex10 -f0 <mat_binaryfile> -ksp_type preonly -pc_type cholesky -mat_type seqsbaij -test_inertia -mat_sigma <sigma>
267:      */
268:     PetscOptionsHasName(PETSC_NULL,"-test_inertia",&flg);
269:     if (flg){
270:       PC        pc;
271:       PetscInt  nneg, nzero, npos;
272:       Mat       F;
273: 
274:       KSPGetPC(ksp,&pc);
275:       PCGetFactoredMatrix(pc,&F);
276:       MatGetInertia(F,&nneg,&nzero,&npos);
277:       PetscPrintf(PETSC_COMM_SELF," MatInertia: nneg: %D, nzero: %D, npos: %D\n",nneg,nzero,npos);
278:     }

280:     /*
281:        Tests "diagonal-scaling of preconditioned residual norm" as used 
282:        by many ODE integrator codes including PVODE. Note this is different
283:        than diagonally scaling the matrix before computing the preconditioner
284:     */
285:     PetscOptionsHasName(PETSC_NULL,"-diagonal_scale",&diagonalscale);
286:     if (diagonalscale) {
287:       PC       pc;
288:       PetscInt j,start,end,n;
289:       Vec      scale;
290: 
291:       KSPGetPC(ksp,&pc);
292:       VecGetSize(x,&n);
293:       VecDuplicate(x,&scale);
294:       VecGetOwnershipRange(scale,&start,&end);
295:       for (j=start; j<end; j++) {
296:         VecSetValue(scale,j,((PetscReal)(j+1))/((PetscReal)n),INSERT_VALUES);
297:       }
298:       VecAssemblyBegin(scale);
299:       VecAssemblyEnd(scale);
300:       PCDiagonalScaleSet(pc,scale);
301:       VecDestroy(scale);

303:     }

305:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
306:                            Solve system
307:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

309:     /*
310:        Begin profiling next stage
311:     */
312:     PreLoadStage("KSPSolve");

314:     /*
315:        Solve linear system; we also explicitly time this stage.
316:     */
317:     PetscGetTime(&tsolve1);
318:     if (trans) {
319:       KSPSolveTranspose(ksp,b,x);
320:       KSPGetIterationNumber(ksp,&its);
321:     } else {
322:       PetscInt  num_rhs=1;
323:       PetscOptionsGetInt(PETSC_NULL,"-num_rhs",&num_rhs,PETSC_NULL);
324:       PetscOptionsHasName(PETSC_NULL,"-cknorm",&cknorm);
325:       while ( num_rhs-- ) {
326:         KSPSolve(ksp,b,x);
327:       }
328:       KSPGetIterationNumber(ksp,&its);
329:       if (cknorm){   /* Check error for each rhs */
330:         if (trans) {
331:           MatMultTranspose(A,x,u);
332:         } else {
333:           MatMult(A,x,u);
334:         }
335:         VecAXPY(&none,b,u);
336:         VecNorm(u,NORM_2,&norm);
337:         PetscPrintf(PETSC_COMM_WORLD,"  Number of iterations = %3D\n",its);
338:         PetscPrintf(PETSC_COMM_WORLD,"  Residual norm %A\n",norm);
339:       }
340:     } /* while ( num_rhs-- ) */
341:     PetscGetTime(&tsolve2);
342:     tsolve = tsolve2 - tsolve1;

344:    /* 
345:        Conclude profiling this stage
346:     */
347:     PreLoadStage("Cleanup");

349:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
350:             Check error, print output, free data structures.
351:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

353:     /* 
354:        Check error
355:     */
356:     if (trans) {
357:       MatMultTranspose(A,x,u);
358:     } else {
359:       MatMult(A,x,u);
360:     }
361:     VecAXPY(&none,b,u);
362:     VecNorm(u,NORM_2,&norm);

364:     /*
365:        Write output (optinally using table for solver details).
366:         - PetscPrintf() handles output for multiprocessor jobs 
367:           by printing from only one processor in the communicator.
368:         - KSPView() prints information about the linear solver.
369:     */
370:     if (table) {
371:       char        *matrixname,kspinfo[120];
372:       PetscViewer viewer;

374:       /*
375:          Open a string viewer; then write info to it.
376:       */
377:       PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
378:       KSPView(ksp,viewer);
379:       PetscStrrchr(file[PreLoadIt],'/',&matrixname);
380:       PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %2.1e %2.1e %2.1e %s \n",
381:                 matrixname,its,norm,tsetup+tsolve,tsetup,tsolve,kspinfo);

383:       /*
384:          Destroy the viewer
385:       */
386:       PetscViewerDestroy(viewer);
387:     } else {
388:       PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
389:       PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A\n",norm);
390:     }

392:     PetscOptionsHasName(PETSC_NULL, "-ksp_reason", &flg);
393:     if (flg){
394:       KSPConvergedReason reason;
395:       KSPGetConvergedReason(ksp,&reason);
396:       PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
397:     }
398: 
399:     } /* while ( num_numfac-- ) */

401:     /* 
402:        Free work space.  All PETSc objects should be destroyed when they
403:        are no longer needed.
404:     */
405:     MatDestroy(A); VecDestroy(b);
406:     VecDestroy(u); VecDestroy(x);
407:     KSPDestroy(ksp);
408:     if (flgB) { MatDestroy(B); }
409:   PreLoadEnd();
410:   /* -----------------------------------------------------------
411:                       End of linear solver loop
412:      ----------------------------------------------------------- */

414:   PetscFinalize();
415:   return 0;
416: }