Actual source code: ex3.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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 10

 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 PetscErrorCode femBg(PetscScalar[][3],PetscScalar*,PetscInt,PetscScalar*,PetscReal);
 40: extern PetscErrorCode femA(AppCtx*,PetscInt,PetscScalar*);
 41: extern PetscErrorCode 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:   PetscOptionsGetInt(NULL,NULL,"-nphase",&nphase,NULL);
 64:   if (nphase > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nphase must be an integer between 1 and 3");

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

 75:   appctx.m          = m;
 76:   appctx.max_probsz = nz;
 77:   appctx.debug      = PETSC_FALSE;
 78:   appctx.useAlhs    = PETSC_FALSE;

 80:   PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);
 81:   PetscOptionsHasName(NULL,NULL,"-useAlhs",&appctx.useAlhs);

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

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

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

103:   /* create the jacobian matrix */
104:   /*----------------------------*/
105:   MatCreate(PETSC_COMM_WORLD, &Jmat);
106:   MatSetSizes(Jmat,PETSC_DECIDE,PETSC_DECIDE,m,m);
107:   MatSetFromOptions(Jmat);
108:   MatSetUp(Jmat);

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

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

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

128:   /* set time-step method */
129:   TSSetType(ts,TSCN);

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

136:   if (appctx.useAlhs) {
137:     /* set the left hand side matrix of Amat*U_t = rhs(U,t) */

139:     /* Note: this approach is incompatible with the finite differenced Jacobian set below because we can't restore the
140:      * Alhs matrix without making a copy.  Either finite difference the entire thing or use analytic Jacobians in both
141:      * places.
142:      */
143:     TSSetIFunction(ts,NULL,TSComputeIFunctionLinear,&appctx);
144:     TSSetIJacobian(ts,appctx.Amat,appctx.Amat,TSComputeIJacobianConstant,&appctx);
145:   }

147:   /* use petsc to compute the jacobian by finite differences */
148:   TSGetSNES(ts,&snes);
149:   SNESSetJacobian(snes,Jmat,Jmat,SNESComputeJacobianDefault,NULL);

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

154: #if defined(PETSC_HAVE_SUNDIALS)
155:   {
156:     TSType    type;
157:     PetscBool sundialstype=PETSC_FALSE;
158:     TSGetType(ts,&type);
159:     PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&sundialstype);
160:     if (sundialstype && appctx.useAlhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use Alhs formulation for TSSUNDIALS type");
161:   }
162: #endif
163:   /* Sets the initial solution */
164:   TSSetSolution(ts,init_sol);

166:   stepsz[0] = 1.0/(2.0*(nz-1)*(nz-1)); /* (mesh_size)^2/2.0 */
167:   ftime     = 0.0;
168:   for (k=0; k<nphase; k++) {
169:     if (nphase > 1) {PetscPrintf(PETSC_COMM_WORLD,"Phase %D initial time %g, stepsz %g, duration: %g\n",k,(double)ftime,(double)stepsz[k],(double)((k+1)*T));}
170:     TSSetInitialTimeStep(ts,ftime,stepsz[k]);
171:     TSSetDuration(ts,max_steps,(k+1)*T);
172:     TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

174:     /* loop over time steps */
175:     /*----------------------*/
176:     TSSolve(ts,init_sol);
177:     TSGetSolveTime(ts,&ftime);
178:     TSGetTimeStepNumber(ts,&steps);
179:     stepsz[k+1] = stepsz[k]*1.5; /* change step size for the next phase */
180:   }

182:   /* free space */
183:   TSDestroy(&ts);
184:   MatDestroy(&appctx.Amat);
185:   MatDestroy(&Jmat);
186:   VecDestroy(&appctx.ksp_rhs);
187:   VecDestroy(&appctx.ksp_sol);
188:   VecDestroy(&init_sol);
189:   VecDestroy(&appctx.solution);
190:   PetscFree(z);

192:   PetscFinalize();
193:   return 0;
194:   }

196: /*------------------------------------------------------------------------
197:   Set exact solution
198:   u(z,t) = sin(6*PI*z)*exp(-36.*PI*PI*t) + 3.*sin(2*PI*z)*exp(-4.*PI*PI*t)
199: --------------------------------------------------------------------------*/
200: PetscScalar exact(PetscScalar z,PetscReal t)
201: {
202:   PetscScalar val, ex1, ex2;

204:   ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t);
205:   ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t);
206:   val = PetscSinScalar(6*PETSC_PI*z)*ex1 + 3.*PetscSinScalar(2*PETSC_PI*z)*ex2;
207:   return val;
208: }

212: /*
213:    Monitor - User-provided routine to monitor the solution computed at
214:    each timestep.  This example plots the solution and computes the
215:    error in two different norms.

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

236:   /* Compute the exact solution */
237:   VecGetArray(appctx->solution,&u_exact);
238:   for (i=0; i<m; i++) u_exact[i] = exact(appctx->z[i+1],time);
239:   VecRestoreArray(appctx->solution,&u_exact);

241:   /* Print debugging information if desired */
242:   if (appctx->debug) {
243:     PetscPrintf(PETSC_COMM_SELF,"Computed solution vector at time %g\n",(double)time);
244:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
245:     PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");
246:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
247:   }

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

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

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

258:   /*
259:      Print debugging information if desired
260:   */
261:   if (appctx->debug) {
262:     PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
263:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
264:   }
265:   return 0;
266: }

268: /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
269: %%      Function to solve a linear system using KSP                                           %%
270: %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/

272: PetscErrorCode Petsc_KSPSolve(AppCtx *obj)
273: {
275:   KSP            ksp;
276:   PC             pc;

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

282:   /*get the preconditioner context, set its type and the tolerances*/
283:   KSPGetPC(ksp,&pc);
284:   PCSetType(pc,PCLU);
285:   KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

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

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

293:   KSPDestroy(&ksp);
294:   return 0;
295: }

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

315: PetscScalar bspl(PetscScalar *x, PetscScalar xx,PetscInt il,PetscInt iq,PetscInt nll[][2],PetscInt id)
316: {
317:   PetscScalar x1,x2,bfcn;
318:   PetscInt    i1,i2,iq1,iq2;

320:   /*** Determine which basis function in interval intrvl is to be used in ***/
321:   iq1 = iq;
322:   if (iq1==0) iq2 = 1;
323:   else iq2 = 0;

325:   /***  Determine endpoint of the interval intrvl ***/
326:   i1=nll[il][iq1];
327:   i2=nll[il][iq2];

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

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

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

361:   /*  quadrature weights  */
362:   qdwt[0] = 1.0/6.0;
363:   qdwt[1] = 4.0/6.0;
364:   qdwt[2] = 1.0/6.0;

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

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

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

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

392:       for (iq=0; iq < 2; iq++) {
393:         ip  = nli[il][iq];
394:         b_z = bspl(z,zz,il,iq,nli,2);
395:         i   = indx[ip];

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

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

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

430:   /*  initializing everything  */

432:   for (i=0; i < nz; i++) {
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++) indx[i]=i-1;
451:   indx[nz-1]=-1;

453:   ipq = 0;

455:   for (il=0; il < nz-1; il++) {
456:     ip           = ipq;
457:     ipq          = ip+1;
458:     zip          = z[ip];
459:     zipq         = z[ipq];
460:     dl           = zipq-zip;
461:     rquad[il][0] = zip;
462:     rquad[il][1] = (0.5)*(zip+zipq);
463:     rquad[il][2] = zipq;
464:     dlen[il]     = fabs(dl);
465:     nli[il][0]   = ip;
466:     nli[il][1]   = ipq;
467:   } /*end for (il)*/

469:   for (il=0; il < nz-1; il++) {
470:     for (iquad=0; iquad < 3; iquad++) {
471:       dd = (dlen[il])*(qdwt[iquad]);
472:       zz = rquad[il][iquad];

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

498: /*---------------------------------------------------------
499:         Function to fill the rhs vector with
500:         By + g values ****
501: ---------------------------------------------------------*/
502: PetscErrorCode rhs(AppCtx *obj,PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t)
503: {
504:   PetscInt       i,j,js,je,jj;
505:   PetscScalar    val,g[num_z],btri[num_z][3],add_term;

508:   for (i=0; i < nz-2; i++) {
509:     for (j=0; j <= 2; j++) btri[i][j]=0.0;
510:     g[i] = 0.0;
511:   }

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

516:   /*  setting the entries of the right hand side vector  */
517:   for (i=0; i < nz-2; i++) {
518:     val = 0.0;
519:     js  = 0;
520:     if (i == 0) js = 1;
521:     je = 2;
522:     if (i == nz-2) je = 1;

524:     for (jj=js; jj <= je; jj++) {
525:       j    = i+jj-1;
526:       val += (btri[i][jj])*(y[j]);
527:     }
528:     add_term = val + g[i];
529:     VecSetValue(obj->ksp_rhs,(PetscInt)i,(PetscScalar)add_term,INSERT_VALUES);
530:   }
531:   VecAssemblyBegin(obj->ksp_rhs);
532:   VecAssemblyEnd(obj->ksp_rhs);

534:   /*  return to main driver function  */
535:   return 0;
536: }

538: /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
539: %%   Function to form the right hand side of the time-stepping problem.                       %%
540: %% -------------------------------------------------------------------------------------------%%
541:   if (useAlhs):
542:     globalout = By+g
543:   else if (!useAlhs):
544:     globalout = f(y,t)=Ainv(By+g),
545:       in which the ksp solver to transform the problem A*ydot=By+g
546:       to the problem ydot=f(y,t)=inv(A)*(By+g)
547: %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/

549: PetscErrorCode RHSfunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
550: {
551:   PetscErrorCode    ierr;
552:   AppCtx            *obj = (AppCtx*)ctx;
553:   PetscScalar       soln[num_z];
554:   const PetscScalar *soln_ptr;
555:   PetscInt          i,nz=obj->nz;
556:   PetscReal         time;

558:   /* get the previous solution to compute updated system */
559:   VecGetArrayRead(globalin,&soln_ptr);
560:   for (i=0; i < num_z-2; i++) soln[i] = soln_ptr[i];
561:   VecRestoreArrayRead(globalin,&soln_ptr);
562:   soln[num_z-1] = 0.0;
563:   soln[num_z-2] = 0.0;

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

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

572:   /* do a ksp solve to get the rhs for the ts problem */
573:   if (obj->useAlhs) {
574:     /* ksp_sol = ksp_rhs */
575:     VecCopy(obj->ksp_rhs,globalout);
576:   } else {
577:     /* ksp_sol = inv(Amat)*ksp_rhs */
578:     Petsc_KSPSolve(obj);
579:     VecCopy(obj->ksp_sol,globalout);
580:   }
581:   return 0;
582: }