Actual source code: ex3.c

petsc-3.3-p7 2013-05-11
  2: static char help[] = "Solves 1D heat equation with FEM formulation.\n\
  3: Input arguments are\n\
  4:   -useAlhs: solve Alhs*U' =  (Arhs*U + g) \n\
  5:             otherwise, solve U' = inv(Alhs)*(Arhs*U + g) \n\n";

  7: /*--------------------------------------------------------------------------
  8:   Solves 1D heat equation U_t = U_xx with FEM formulation:
  9:                           Alhs*U' = rhs (= Arhs*U + g)
 10:   We thank Chris Cox <clcox@clemson.edu> for contributing the original code
 11: ----------------------------------------------------------------------------*/

 13: #include <petscksp.h>
 14: #include <petscts.h>

 16: /* special variable - max size of all arrays  */
 17: #define num_z 60

 19: /* 
 20:    User-defined application context - contains data needed by the 
 21:    application-provided call-back routines.
 22: */
 23: typedef struct{
 24:   Mat            Amat;            /* left hand side matrix */
 25:   Vec            ksp_rhs,ksp_sol; /* working vectors for formulating inv(Alhs)*(Arhs*U+g) */
 26:   int            max_probsz;      /* max size of the problem */
 27:   PetscBool      useAlhs;         /* flag (1 indicates solving Alhs*U' = Arhs*U+g */
 28:   int            nz;              /* total number of grid points */
 29:   PetscInt       m;               /* total number of interio grid points */
 30:   Vec            solution;        /* global exact ts solution vector */
 31:   PetscScalar    *z;              /* array of grid points */
 32:   PetscBool      debug;           /* flag (1 indicates activation of debugging printouts) */
 33: } AppCtx;

 35: extern PetscScalar exact(PetscScalar,PetscReal);
 36: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 37: extern PetscErrorCode Petsc_KSPSolve(AppCtx*);
 38: extern PetscScalar bspl(PetscScalar*,PetscScalar,PetscInt,PetscInt,PetscInt[][2],PetscInt);
 39: extern void femBg(PetscScalar[][3],PetscScalar*,PetscInt,PetscScalar*,PetscReal);
 40: extern void femA(AppCtx*,PetscInt,PetscScalar*);
 41: extern void rhs(AppCtx*,PetscScalar*, PetscInt, PetscScalar*,PetscReal);
 42: extern PetscErrorCode RHSfunction(TS,PetscReal,Vec,Vec,void*);

 46: int main(int argc,char **argv)
 47: {
 48:   PetscInt       i,m,nz,steps,max_steps,k,nphase=1;
 49:   PetscScalar    zInitial,zFinal,val,*z;
 50:   PetscReal      stepsz[4],T,ftime;
 52:   TS             ts;
 53:   SNES           snes;
 54:   Mat            Jmat;
 55:   AppCtx         appctx;   /* user-defined application context */
 56:   Vec            init_sol; /* ts solution vector */
 57:   PetscMPIInt    size;

 59:   PetscInitialize(&argc,&argv,(char*)0,help);
 60:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 61:   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only");

 63:   PetscOptionsHasName(PETSC_NULL,"-debug",&appctx.debug);
 64:   PetscOptionsHasName(PETSC_NULL,"-useAlhs",&appctx.useAlhs);
 65:   PetscOptionsGetInt(PETSC_NULL,"-nphase",&nphase,PETSC_NULL);
 66:   if (nphase > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nphase must be an integer between 1 and 3");

 68:   /* initializations */
 69:   zInitial       = 0.0;
 70:   zFinal         = 1.0;
 71:   T              = 0.014/nphase;
 72:   nz             = num_z;
 73:   m              = nz-2;
 74:   appctx.nz      = nz;
 75:   max_steps      = (PetscInt)10000;

 77:   appctx.m          = m;
 78:   appctx.max_probsz = nz;
 79:   appctx.debug      = PETSC_FALSE;
 80:   appctx.useAlhs    = PETSC_FALSE;

 82:   /* create vector to hold ts solution */
 83:   /*-----------------------------------*/
 84:   VecCreate(PETSC_COMM_WORLD, &init_sol);
 85:   VecSetSizes(init_sol, PETSC_DECIDE, m);
 86:   VecSetFromOptions(init_sol);

 88:   /* create vector to hold true ts soln for comparison */
 89:   VecDuplicate(init_sol, &appctx.solution);

 91:   /* create LHS matrix Amat */
 92:   /*------------------------*/
 93:   MatCreateSeqAIJ(PETSC_COMM_WORLD, m, m, 3, PETSC_NULL, &appctx.Amat);
 94:   MatSetFromOptions(appctx.Amat);
 95:   /* set space grid points - interio points only! */
 96:   PetscMalloc((nz+1)*sizeof(PetscScalar),&z);
 97:   for (i=0; i<nz; i++) z[i]=(i)*((zFinal-zInitial)/(nz-1));
 98:   appctx.z = z;
 99:   femA(&appctx,nz,z);

101:   /* create the jacobian matrix */
102:   /*----------------------------*/
103:   MatCreate(PETSC_COMM_WORLD, &Jmat);
104:   MatSetSizes(Jmat,PETSC_DECIDE,PETSC_DECIDE,m,m);
105:   MatSetFromOptions(Jmat);

107:   /* create working vectors for formulating rhs=inv(Alhs)*(Arhs*U + g) */
108:   VecDuplicate(init_sol,&appctx.ksp_rhs);
109:   VecDuplicate(init_sol,&appctx.ksp_sol);

111:   /* set intial guess */
112:   /*------------------*/
113:   for(i=0; i<nz-2; i++){
114:     val = exact(z[i+1], 0.0);
115:     VecSetValue(init_sol,i,(PetscScalar)val,INSERT_VALUES);
116:   }
117:   VecAssemblyBegin(init_sol);
118:   VecAssemblyEnd(init_sol);

120:   /*create a time-stepping context and set the problem type */
121:   /*--------------------------------------------------------*/
122:   TSCreate(PETSC_COMM_WORLD, &ts);
123:   TSSetProblemType(ts,TS_NONLINEAR);

125:   /* set time-step method */
126:   TSSetType(ts,TSCN);

128:   /* Set optional user-defined monitoring routine */
129:   TSMonitorSet(ts,Monitor,&appctx,PETSC_NULL);
130:   /* set the right hand side of U_t = RHSfunction(U,t) */
131:   TSSetRHSFunction(ts,PETSC_NULL,(PetscErrorCode (*)(TS,PetscScalar,Vec,Vec,void*))RHSfunction,&appctx);

133:   if (appctx.useAlhs){
134:     /* set the left hand side matrix of Amat*U_t = rhs(U,t) */
135:     TSSetIFunction(ts,PETSC_NULL,TSComputeIFunctionLinear,&appctx);
136:     TSSetIJacobian(ts,appctx.Amat,appctx.Amat,TSComputeIJacobianConstant,&appctx);
137:   }

139:   /* use petsc to compute the jacobian by finite differences */
140:   TSGetSNES(ts,&snes);
141:   SNESSetJacobian(snes,Jmat,Jmat,SNESDefaultComputeJacobian,PETSC_NULL);

143:   /* get the command line options if there are any and set them */
144:   TSSetFromOptions(ts);

146: #ifdef PETSC_HAVE_SUNDIALS
147:   {
148:     const TSType   type;
149:     PetscBool      sundialstype=PETSC_FALSE;
150:     TSGetType(ts,&type);
151:     PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&sundialstype);
152:     if (sundialstype && appctx.useAlhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use Alhs formulation for TSSUNDIALS type");
153:   }
154: #endif
155:   /* Sets the initial solution */
156:   TSSetSolution(ts,init_sol);

158:   stepsz[0] = 1.0/(2.0*(nz-1)*(nz-1)); /* (mesh_size)^2/2.0 */
159:   ftime = 0.0;
160:   for (k=0; k<nphase; k++){
161:     if (nphase > 1) {
162:       printf("Phase %d: initial time %g, stepsz %g, duration: %g\n",k,ftime,stepsz[k],(k+1)*T);
163:     }
164:     TSSetInitialTimeStep(ts,ftime,stepsz[k]);
165:     TSSetDuration(ts,max_steps,(k+1)*T);

167:     /* loop over time steps */
168:     /*----------------------*/
169:     TSSolve(ts,init_sol,&ftime);
170:     TSGetTimeStepNumber(ts,&steps);
171:     stepsz[k+1] = stepsz[k]*1.5; /* change step size for the next phase */
172:   }

174:   /* free space */
175:   TSDestroy(&ts);
176:   MatDestroy(&appctx.Amat);
177:   MatDestroy(&Jmat);
178:   VecDestroy(&appctx.ksp_rhs);
179:   VecDestroy(&appctx.ksp_sol);
180:   VecDestroy(&init_sol);
181:   VecDestroy(&appctx.solution);
182:   PetscFree(z);

184:   PetscFinalize();
185:   return 0;
186: }

188: /*------------------------------------------------------------------------
189:   Set exact solution 
190:   u(z,t) = sin(6*PI*z)*exp(-36.*PI*PI*t) + 3.*sin(2*PI*z)*exp(-4.*PI*PI*t)
191: --------------------------------------------------------------------------*/
192: PetscScalar exact(PetscScalar z,PetscReal t)
193: {
194:   PetscScalar val, ex1, ex2;

196:   ex1 = exp(-36.*PETSC_PI*PETSC_PI*t);
197:   ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
198:   val = sin(6*PETSC_PI*z)*ex1 + 3.*sin(2*PETSC_PI*z)*ex2;
199:   return val;
200: }

204: /*
205:    Monitor - User-provided routine to monitor the solution computed at 
206:    each timestep.  This example plots the solution and computes the
207:    error in two different norms.

209:    Input Parameters:
210:    ts     - the timestep context
211:    step   - the count of the current step (with 0 meaning the
212:              initial condition)
213:    time   - the current time
214:    u      - the solution at this timestep
215:    ctx    - the user-provided context for this monitoring routine.
216:             In this case we use the application context which contains 
217:             information about the problem size, workspace and the exact 
218:             solution.
219: */
220: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
221: {
222:   AppCtx         *appctx = (AppCtx*)ctx;
224:   PetscInt       i,m=appctx->m;
225:   PetscReal      norm_2,norm_max,h=1.0/(m+1);
226:   PetscScalar    *u_exact;

228:   /* Compute the exact solution */
229:   VecGetArray(appctx->solution,&u_exact);
230:   for (i=0; i<m; i++){
231:     u_exact[i] = exact(appctx->z[i+1],time);
232:   }
233:   VecRestoreArray(appctx->solution,&u_exact);

235:   /* Print debugging information if desired */
236:   if (appctx->debug) {
237:      PetscPrintf(PETSC_COMM_SELF,"Computed solution vector at time %g\n",time);
238:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
239:      PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");
240:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
241:   }

243:   /* Compute the 2-norm and max-norm of the error */
244:   VecAXPY(appctx->solution,-1.0,u);
245:   VecNorm(appctx->solution,NORM_2,&norm_2);

247:   norm_2 = PetscSqrtReal(h)*norm_2;
248:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

250:   PetscPrintf(PETSC_COMM_SELF,"Timestep %D: time = %G, 2-norm error = %6.4f, max norm error = %6.4f\n",
251:               step,time,norm_2,norm_max);

253:   /*
254:      Print debugging information if desired
255:   */
256:   if (appctx->debug) {
257:      PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
258:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
259:   }
260:   return 0;
261: }

263: /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
264: %%      Function to solve a linear system using KSP                                           %%
265: %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/

267: PetscErrorCode Petsc_KSPSolve(AppCtx *obj)
268: {
269:    PetscErrorCode  ierr;
270:    KSP             ksp;
271:    PC              pc;

273:    /*create the ksp context and set the operators,that is, associate the system matrix with it*/
274:    KSPCreate(PETSC_COMM_WORLD,&ksp);
275:    KSPSetOperators(ksp,obj->Amat,obj->Amat,DIFFERENT_NONZERO_PATTERN);

277:    /*get the preconditioner context, set its type and the tolerances*/
278:    KSPGetPC(ksp,&pc);
279:    PCSetType(pc,PCLU);
280:    KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

282:    /*get the command line options if there are any and set them*/
283:    KSPSetFromOptions(ksp);

285:    /*get the linear system (ksp) solve*/
286:    KSPSolve(ksp,obj->ksp_rhs,obj->ksp_sol);

288:    KSPDestroy(&ksp);
289:    return 0;
290: }

292: /***********************************************************************
293:  * Function to return value of basis function or derivative of basis   *
294:  *              function.                                              *
295:  ***********************************************************************
296:  *                                                                     *
297:  *       Arguments:                                                    *
298:  *         x       = array of xpoints or nodal values                  *
299:  *         xx      = point at which the basis function is to be        *
300:  *                     evaluated.                                      *
301:  *         il      = interval containing xx.                           *
302:  *         iq      = indicates which of the two basis functions in     *
303:  *                     interval intrvl should be used                  *
304:  *         nll     = array containing the endpoints of each interval.  *
305:  *         id      = If id ~= 2, the value of the basis function       *
306:  *                     is calculated; if id = 2, the value of the      *
307:  *                     derivative of the basis function is returned.   *
308:  ***********************************************************************/

310: PetscScalar bspl(PetscScalar *x, PetscScalar xx,PetscInt il,PetscInt iq,PetscInt nll[][2],PetscInt id)
311: {
312:   PetscScalar  x1,x2,bfcn;
313:   PetscInt     i1,i2,iq1,iq2;

315:   /*** Determine which basis function in interval intrvl is to be used in ***/
316:   iq1 = iq;
317:   if(iq1==0) iq2 = 1;
318:   else iq2 = 0;

320:   /***  Determine endpoint of the interval intrvl ***/
321:   i1=nll[il][iq1];
322:   i2=nll[il][iq2];

324:   /*** Determine nodal values at the endpoints of the interval intrvl ***/
325:   x1=x[i1];
326:   x2=x[i2];
327:   //printf("x1=%g\tx2=%g\txx=%g\n",x1,x2,xx);
328:   /*** Evaluate basis function ***/
329:   if(id == 2) bfcn=(1.0)/(x1-x2);
330:   else bfcn=(xx-x2)/(x1-x2);
331:   //printf("bfcn=%g\n",bfcn);
332:   return bfcn;
333: }

335: /*---------------------------------------------------------
336:   Function called by rhs function to get B and g 
337: ---------------------------------------------------------*/
338: void femBg(PetscScalar btri[][3],PetscScalar *f,PetscInt nz,PetscScalar *z, PetscReal t)
339: {
340:   PetscInt    i,j,jj,il,ip,ipp,ipq,iq,iquad,iqq;
341:   PetscInt    nli[num_z][2],indx[num_z];
342:   PetscScalar dd,dl,zip,zipq,zz,bb,b_z,bbb,bb_z,bij;
343:   PetscScalar zquad[num_z][3],dlen[num_z],qdwt[3];

345:   /*  initializing everything - btri and f are initialized in rhs.c  */
346:   for(i=0; i < nz; i++){
347:     nli[i][0] = 0;
348:     nli[i][1] = 0;
349:     indx[i] = 0;
350:     zquad[i][0] = 0.0;
351:     zquad[i][1] = 0.0;
352:     zquad[i][2] = 0.0;
353:     dlen[i] = 0.0;
354:   }/*end for(i)*/

356:   /*  quadrature weights  */
357:   qdwt[0] = 1.0/6.0;
358:   qdwt[1] = 4.0/6.0;
359:   qdwt[2] = 1.0/6.0;

361:   /* 1st and last nodes have Dirichlet boundary condition -
362:      set indices there to -1 */

364:   for(i=0; i < nz-1; i++){
365:     indx[i]=i-1;
366:   }
367:   indx[nz-1]=-1;

369:   ipq = 0;
370:   for (il=0; il < nz-1; il++){
371:     ip = ipq;
372:     ipq = ip+1;
373:     zip = z[ip];
374:     zipq = z[ipq];
375:     dl = zipq-zip;
376:     zquad[il][0] = zip;
377:     zquad[il][1] = (0.5)*(zip+zipq);
378:     zquad[il][2] = zipq;
379:     dlen[il] = fabs(dl);
380:     nli[il][0] = ip;
381:     nli[il][1] = ipq;
382:   }

384:   for (il=0; il < nz-1; il++){
385:     for (iquad=0; iquad < 3; iquad++){
386:       dd = (dlen[il])*(qdwt[iquad]);
387:       zz = zquad[il][iquad];

389:       for (iq=0; iq < 2; iq++){
390:         ip = nli[il][iq];
391:         bb = bspl(z,zz,il,iq,nli,1);
392:         b_z = bspl(z,zz,il,iq,nli,2);
393:         i = indx[ip];

395:         if(i > -1){
396:           for(iqq=0; iqq < 2; iqq++){
397:             ipp = nli[il][iqq];
398:             bbb = bspl(z,zz,il,iqq,nli,1);
399:             bb_z = bspl(z,zz,il,iqq,nli,2);
400:             j = indx[ipp];
401:             bij = -b_z*bb_z;

403:             if (j > -1){
404:               jj = 1+j-i;
405:               btri[i][jj] += bij*dd;
406:             } else {
407:               f[i] += bij*dd*exact(z[ipp], t);
408:               // f[i] += 0.0;
409:               // if(il==0 && j==-1){
410:               // f[i] += bij*dd*exact(zz,t);
411:               // }/*end if*/
412:             } /*end else*/
413:           }/*end for(iqq)*/
414:         }/*end if(i>0)*/
415:       }/*end for(iq)*/
416:     }/*end for(iquad)*/
417:   }/*end for(il)*/
418:   return;
419: }

421: void femA(AppCtx *obj,PetscInt nz,PetscScalar *z)
422: {
423:   PetscInt        i,j,il,ip,ipp,ipq,iq,iquad,iqq;
424:   PetscInt        nli[num_z][2],indx[num_z];
425:   PetscScalar     dd,dl,zip,zipq,zz,bb,bbb,aij;
426:   PetscScalar     rquad[num_z][3],dlen[num_z],qdwt[3],add_term;
427:   PetscErrorCode  ierr;

429:   /*  initializing everything  */

431:   for(i=0; i < nz; i++)
432:   {
433:      nli[i][0] = 0;
434:      nli[i][1] = 0;
435:      indx[i] = 0;
436:      rquad[i][0] = 0.0;
437:      rquad[i][1] = 0.0;
438:      rquad[i][2] = 0.0;
439:      dlen[i] = 0.0;
440:   }/*end for(i)*/

442:   /*  quadrature weights  */
443:   qdwt[0] = 1.0/6.0;
444:   qdwt[1] = 4.0/6.0;
445:   qdwt[2] = 1.0/6.0;

447:   /* 1st and last nodes have Dirichlet boundary condition - 
448:      set indices there to -1 */

450:   for(i=0; i < nz-1; i++)
451:   {
452:      indx[i]=i-1;

454:   }/*end for(i)*/
455:   indx[nz-1]=-1;

457:   ipq = 0;

459:   for(il=0; il < nz-1; il++)
460:   {
461:      ip = ipq;
462:      ipq = ip+1;
463:      zip = z[ip];
464:      zipq = z[ipq];
465:      dl = zipq-zip;
466:      rquad[il][0] = zip;
467:      rquad[il][1] = (0.5)*(zip+zipq);
468:      rquad[il][2] = zipq;
469:      dlen[il] = fabs(dl);
470:      nli[il][0] = ip;
471:      nli[il][1] = ipq;

473:   }/*end for(il)*/

475:   for(il=0; il < nz-1; il++){
476:     for(iquad=0; iquad < 3; iquad++){
477:       dd = (dlen[il])*(qdwt[iquad]);
478:       zz = rquad[il][iquad];

480:       for(iq=0; iq < 2; iq++){
481:         ip = nli[il][iq];
482:         bb = bspl(z,zz,il,iq,nli,1);
483:         i = indx[ip];
484:         if(i > -1){
485:           for(iqq=0; iqq < 2; iqq++){
486:             ipp = nli[il][iqq];
487:             bbb = bspl(z,zz,il,iqq,nli,1);
488:             j = indx[ipp];
489:             aij = bb*bbb;
490:             if(j > -1) {
491:               add_term = aij*dd;
492:               MatSetValue(obj->Amat,i,j,add_term,ADD_VALUES);
493:             }/*endif*/
494:           }/*end for(iqq)*/
495:         }/*end if(i>0)*/
496:       }/*end for(iq)*/
497:     }/*end for(iquad)*/
498:   }/*end for(il)*/
499:   MatAssemblyBegin(obj->Amat,MAT_FINAL_ASSEMBLY);
500:   MatAssemblyEnd(obj->Amat,MAT_FINAL_ASSEMBLY);
501:   return;
502: }

504: /*---------------------------------------------------------
505:         Function to fill the rhs vector with
506:         By + g values ****
507: ---------------------------------------------------------*/
508: void rhs(AppCtx *obj,PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t)
509: {
510:   PetscInt        i,j,js,je,jj;
511:   PetscScalar     val,g[num_z],btri[num_z][3],add_term;
512:   PetscErrorCode  ierr;

514:   for(i=0; i < nz-2; i++){
515:     for(j=0; j <= 2; j++){
516:       btri[i][j]=0.0;
517:     }
518:     g[i] = 0.0;
519:   }

521:   /*  call femBg to set the tri-diagonal b matrix and vector g  */
522:   femBg(btri,g,nz,z,t);

524:   /*  setting the entries of the right hand side vector  */
525:   for(i=0; i < nz-2; i++){
526:     val = 0.0;
527:     js = 0;
528:     if(i == 0) js = 1;
529:     je = 2;
530:     if(i == nz-2) je = 1;

532:     for(jj=js; jj <= je; jj++){
533:       j = i+jj-1;
534:       val += (btri[i][jj])*(y[j]);
535:     }
536:     add_term = val + g[i];
537:     VecSetValue(obj->ksp_rhs,(PetscInt)i,(PetscScalar)add_term,INSERT_VALUES);
538:   }
539:   VecAssemblyBegin(obj->ksp_rhs);
540:   VecAssemblyEnd(obj->ksp_rhs);

542:   /*  return to main driver function  */
543:   return;
544: }

546: /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
547: %%   Function to form the right hand side of the time-stepping problem.                       %%
548: %% -------------------------------------------------------------------------------------------%%
549:   if (useAlhs):
550:     globalout = By+g
551:   else if (!useAlhs):
552:     globalout = f(y,t)=Ainv(By+g),
553:       in which the ksp solver to transform the problem A*ydot=By+g
554:       to the problem ydot=f(y,t)=inv(A)*(By+g)
555: %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/

557: PetscErrorCode RHSfunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
558: {
560:   AppCtx         *obj = (AppCtx*)ctx;
561:   PetscScalar    *soln_ptr,soln[num_z-2];
562:   PetscInt       i,nz=obj->nz;
563:   PetscReal      time;

565:   /* get the previous solution to compute updated system */
566:   VecGetArray(globalin,&soln_ptr);
567:   for(i=0;i < num_z-2;i++){
568:     soln[i] = soln_ptr[i];
569:   }
570:   VecRestoreArray(globalin,&soln_ptr);

572:   /* clear out the matrix and rhs for ksp to keep things straight */
573:   VecSet(obj->ksp_rhs,(PetscScalar)0.0);

575:   time = t;
576:   /* get the updated system */
577:   rhs(obj,soln,nz,obj->z,time); /* setup of the By+g rhs */

579:   /* do a ksp solve to get the rhs for the ts problem */
580:   if (obj->useAlhs){
581:     /* ksp_sol = ksp_rhs */
582:     VecCopy(obj->ksp_rhs,globalout);
583:   } else {
584:     /* ksp_sol = inv(Amat)*ksp_rhs */
585:     Petsc_KSPSolve(obj);
586:     VecCopy(obj->ksp_sol,globalout);
587:   }
588:   return 0;
589: }