Actual source code: ex3.c

petsc-3.4.2 2013-07-02
  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(NULL,"-debug",&appctx.debug);
 64:   PetscOptionsHasName(NULL,"-useAlhs",&appctx.useAlhs);
 65:   PetscOptionsGetInt(NULL,"-nphase",&nphase,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, 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,NULL);
130:   /* set the right hand side of U_t = RHSfunction(U,t) */
131:   TSSetRHSFunction(ts,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,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,SNESComputeJacobianDefault,NULL);

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

146: #if defined(PETSC_HAVE_SUNDIALS)
147:   {
148:     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) printf("Phase %d: initial time %g, stepsz %g, duration: %g\n",k,ftime,stepsz[k],(k+1)*T);
162:     TSSetInitialTimeStep(ts,ftime,stepsz[k]);
163:     TSSetDuration(ts,max_steps,(k+1)*T);

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

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

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

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

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

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

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

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

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

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

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

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

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

260: /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
261: %%      Function to solve a linear system using KSP                                           %%
262: %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/

264: PetscErrorCode Petsc_KSPSolve(AppCtx *obj)
265: {
267:   KSP            ksp;
268:   PC             pc;

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

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

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

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

285:   KSPDestroy(&ksp);
286:   return 0;
287: }

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

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

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

317:   /***  Determine endpoint of the interval intrvl ***/
318:   i1=nll[il][iq1];
319:   i2=nll[il][iq2];

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

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

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

353:   /*  quadrature weights  */
354:   qdwt[0] = 1.0/6.0;
355:   qdwt[1] = 4.0/6.0;
356:   qdwt[2] = 1.0/6.0;

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

361:   for (i=0; i < nz-1; i++) indx[i] = i-1;
362:   indx[nz-1] = -1;

364:   ipq = 0;
365:   for (il=0; il < nz-1; il++) {
366:     ip           = ipq;
367:     ipq          = ip+1;
368:     zip          = z[ip];
369:     zipq         = z[ipq];
370:     dl           = zipq-zip;
371:     zquad[il][0] = zip;
372:     zquad[il][1] = (0.5)*(zip+zipq);
373:     zquad[il][2] = zipq;
374:     dlen[il]     = fabs(dl);
375:     nli[il][0]   = ip;
376:     nli[il][1]   = ipq;
377:   }

379:   for (il=0; il < nz-1; il++) {
380:     for (iquad=0; iquad < 3; iquad++) {
381:       dd = (dlen[il])*(qdwt[iquad]);
382:       zz = zquad[il][iquad];

384:       for (iq=0; iq < 2; iq++) {
385:         ip  = nli[il][iq];
386:         bb  = bspl(z,zz,il,iq,nli,1);
387:         b_z = bspl(z,zz,il,iq,nli,2);
388:         i   = indx[ip];

390:         if (i > -1) {
391:           for (iqq=0; iqq < 2; iqq++) {
392:             ipp  = nli[il][iqq];
393:             bbb  = bspl(z,zz,il,iqq,nli,1);
394:             bb_z = bspl(z,zz,il,iqq,nli,2);
395:             j    = indx[ipp];
396:             bij  = -b_z*bb_z;

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

416: void femA(AppCtx *obj,PetscInt nz,PetscScalar *z)
417: {
418:   PetscInt       i,j,il,ip,ipp,ipq,iq,iquad,iqq;
419:   PetscInt       nli[num_z][2],indx[num_z];
420:   PetscScalar    dd,dl,zip,zipq,zz,bb,bbb,aij;
421:   PetscScalar    rquad[num_z][3],dlen[num_z],qdwt[3],add_term;

424:   /*  initializing everything  */

426:   for (i=0; i < nz; i++) {
427:     nli[i][0]   = 0;
428:     nli[i][1]   = 0;
429:     indx[i]     = 0;
430:     rquad[i][0] = 0.0;
431:     rquad[i][1] = 0.0;
432:     rquad[i][2] = 0.0;
433:     dlen[i]     = 0.0;
434:   } /*end for (i)*/

436:   /*  quadrature weights  */
437:   qdwt[0] = 1.0/6.0;
438:   qdwt[1] = 4.0/6.0;
439:   qdwt[2] = 1.0/6.0;

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

444:   for (i=0; i < nz-1; i++) indx[i]=i-1;
445:   indx[nz-1]=-1;

447:   ipq = 0;

449:   for (il=0; il < nz-1; il++) {
450:     ip           = ipq;
451:     ipq          = ip+1;
452:     zip          = z[ip];
453:     zipq         = z[ipq];
454:     dl           = zipq-zip;
455:     rquad[il][0] = zip;
456:     rquad[il][1] = (0.5)*(zip+zipq);
457:     rquad[il][2] = zipq;
458:     dlen[il]     = fabs(dl);
459:     nli[il][0]   = ip;
460:     nli[il][1]   = ipq;
461:   } /*end for (il)*/

463:   for (il=0; il < nz-1; il++) {
464:     for (iquad=0; iquad < 3; iquad++) {
465:       dd = (dlen[il])*(qdwt[iquad]);
466:       zz = rquad[il][iquad];

468:       for (iq=0; iq < 2; iq++) {
469:         ip = nli[il][iq];
470:         bb = bspl(z,zz,il,iq,nli,1);
471:         i = indx[ip];
472:         if (i > -1) {
473:           for (iqq=0; iqq < 2; iqq++) {
474:             ipp = nli[il][iqq];
475:             bbb = bspl(z,zz,il,iqq,nli,1);
476:             j = indx[ipp];
477:             aij = bb*bbb;
478:             if (j > -1) {
479:               add_term = aij*dd;
480:               MatSetValue(obj->Amat,i,j,add_term,ADD_VALUES);
481:             }/*endif*/
482:           } /*end for (iqq)*/
483:         } /*end if (i>0)*/
484:       } /*end for (iq)*/
485:     } /*end for (iquad)*/
486:   } /*end for (il)*/
487:   MatAssemblyBegin(obj->Amat,MAT_FINAL_ASSEMBLY);
488:   MatAssemblyEnd(obj->Amat,MAT_FINAL_ASSEMBLY);
489:   return;
490: }

492: /*---------------------------------------------------------
493:         Function to fill the rhs vector with
494:         By + g values ****
495: ---------------------------------------------------------*/
496: void rhs(AppCtx *obj,PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t)
497: {
498:   PetscInt       i,j,js,je,jj;
499:   PetscScalar    val,g[num_z],btri[num_z][3],add_term;

502:   for (i=0; i < nz-2; i++) {
503:     for (j=0; j <= 2; j++) btri[i][j]=0.0;
504:     g[i] = 0.0;
505:   }

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

510:   /*  setting the entries of the right hand side vector  */
511:   for (i=0; i < nz-2; i++) {
512:     val = 0.0;
513:     js  = 0;
514:     if (i == 0) js = 1;
515:     je = 2;
516:     if (i == nz-2) je = 1;

518:     for (jj=js; jj <= je; jj++) {
519:       j    = i+jj-1;
520:       val += (btri[i][jj])*(y[j]);
521:     }
522:     add_term = val + g[i];
523:     VecSetValue(obj->ksp_rhs,(PetscInt)i,(PetscScalar)add_term,INSERT_VALUES);
524:   }
525:   VecAssemblyBegin(obj->ksp_rhs);
526:   VecAssemblyEnd(obj->ksp_rhs);

528:   /*  return to main driver function  */
529:   return;
530: }

532: /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
533: %%   Function to form the right hand side of the time-stepping problem.                       %%
534: %% -------------------------------------------------------------------------------------------%%
535:   if (useAlhs):
536:     globalout = By+g
537:   else if (!useAlhs):
538:     globalout = f(y,t)=Ainv(By+g),
539:       in which the ksp solver to transform the problem A*ydot=By+g
540:       to the problem ydot=f(y,t)=inv(A)*(By+g)
541: %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/

543: PetscErrorCode RHSfunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
544: {
546:   AppCtx         *obj = (AppCtx*)ctx;
547:   PetscScalar    *soln_ptr,soln[num_z-2];
548:   PetscInt       i,nz=obj->nz;
549:   PetscReal      time;

551:   /* get the previous solution to compute updated system */
552:   VecGetArray(globalin,&soln_ptr);
553:   for (i=0; i < num_z-2; i++) soln[i] = soln_ptr[i];
554:   VecRestoreArray(globalin,&soln_ptr);

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

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

563:   /* do a ksp solve to get the rhs for the ts problem */
564:   if (obj->useAlhs) {
565:     /* ksp_sol = ksp_rhs */
566:     VecCopy(obj->ksp_rhs,globalout);
567:   } else {
568:     /* ksp_sol = inv(Amat)*ksp_rhs */
569:     Petsc_KSPSolve(obj);
570:     VecCopy(obj->ksp_sol,globalout);
571:   }
572:   return 0;
573: }