xref: /petsc/src/ts/tests/ex3.c (revision d2522c19e8fa9bca20aaca277941d9a63e71db6a)
1 
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";
6 
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 ----------------------------------------------------------------------------*/
12 
13 #include <petscksp.h>
14 #include <petscts.h>
15 
16 /* special variable - max size of all arrays  */
17 #define num_z 10
18 
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;
34 
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 *);
43 
44 int main(int argc, char **argv) {
45   PetscInt    i, m, nz, steps, max_steps, k, nphase = 1;
46   PetscScalar zInitial, zFinal, val, *z;
47   PetscReal   stepsz[4], T, ftime;
48   TS          ts;
49   SNES        snes;
50   Mat         Jmat;
51   AppCtx      appctx;   /* user-defined application context */
52   Vec         init_sol; /* ts solution vector */
53   PetscMPIInt size;
54 
55   PetscFunctionBeginUser;
56   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
57   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
58   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
59 
60   /* initializations */
61   zInitial  = 0.0;
62   zFinal    = 1.0;
63   nz        = num_z;
64   m         = nz - 2;
65   appctx.nz = nz;
66   max_steps = (PetscInt)10000;
67 
68   appctx.m          = m;
69   appctx.max_probsz = nz;
70   appctx.debug      = PETSC_FALSE;
71   appctx.useAlhs    = PETSC_FALSE;
72 
73   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "", "");
74   PetscCall(PetscOptionsName("-debug", NULL, NULL, &appctx.debug));
75   PetscCall(PetscOptionsName("-useAlhs", NULL, NULL, &appctx.useAlhs));
76   PetscCall(PetscOptionsRangeInt("-nphase", NULL, NULL, nphase, &nphase, NULL, 1, 3));
77   PetscOptionsEnd();
78   T = 0.014 / nphase;
79 
80   /* create vector to hold ts solution */
81   /*-----------------------------------*/
82   PetscCall(VecCreate(PETSC_COMM_WORLD, &init_sol));
83   PetscCall(VecSetSizes(init_sol, PETSC_DECIDE, m));
84   PetscCall(VecSetFromOptions(init_sol));
85 
86   /* create vector to hold true ts soln for comparison */
87   PetscCall(VecDuplicate(init_sol, &appctx.solution));
88 
89   /* create LHS matrix Amat */
90   /*------------------------*/
91   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, m, m, 3, NULL, &appctx.Amat));
92   PetscCall(MatSetFromOptions(appctx.Amat));
93   PetscCall(MatSetUp(appctx.Amat));
94   /* set space grid points - interio points only! */
95   PetscCall(PetscMalloc1(nz + 1, &z));
96   for (i = 0; i < nz; i++) z[i] = (i) * ((zFinal - zInitial) / (nz - 1));
97   appctx.z = z;
98   femA(&appctx, nz, z);
99 
100   /* create the jacobian matrix */
101   /*----------------------------*/
102   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jmat));
103   PetscCall(MatSetSizes(Jmat, PETSC_DECIDE, PETSC_DECIDE, m, m));
104   PetscCall(MatSetFromOptions(Jmat));
105   PetscCall(MatSetUp(Jmat));
106 
107   /* create working vectors for formulating rhs=inv(Alhs)*(Arhs*U + g) */
108   PetscCall(VecDuplicate(init_sol, &appctx.ksp_rhs));
109   PetscCall(VecDuplicate(init_sol, &appctx.ksp_sol));
110 
111   /* set initial guess */
112   /*-------------------*/
113   for (i = 0; i < nz - 2; i++) {
114     val = exact(z[i + 1], 0.0);
115     PetscCall(VecSetValue(init_sol, i, (PetscScalar)val, INSERT_VALUES));
116   }
117   PetscCall(VecAssemblyBegin(init_sol));
118   PetscCall(VecAssemblyEnd(init_sol));
119 
120   /*create a time-stepping context and set the problem type */
121   /*--------------------------------------------------------*/
122   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
123   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
124 
125   /* set time-step method */
126   PetscCall(TSSetType(ts, TSCN));
127 
128   /* Set optional user-defined monitoring routine */
129   PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
130   /* set the right hand side of U_t = RHSfunction(U,t) */
131   PetscCall(TSSetRHSFunction(ts, NULL, (PetscErrorCode(*)(TS, PetscScalar, Vec, Vec, void *))RHSfunction, &appctx));
132 
133   if (appctx.useAlhs) {
134     /* set the left hand side matrix of Amat*U_t = rhs(U,t) */
135 
136     /* Note: this approach is incompatible with the finite differenced Jacobian set below because we can't restore the
137      * Alhs matrix without making a copy.  Either finite difference the entire thing or use analytic Jacobians in both
138      * places.
139      */
140     PetscCall(TSSetIFunction(ts, NULL, TSComputeIFunctionLinear, &appctx));
141     PetscCall(TSSetIJacobian(ts, appctx.Amat, appctx.Amat, TSComputeIJacobianConstant, &appctx));
142   }
143 
144   /* use petsc to compute the jacobian by finite differences */
145   PetscCall(TSGetSNES(ts, &snes));
146   PetscCall(SNESSetJacobian(snes, Jmat, Jmat, SNESComputeJacobianDefault, NULL));
147 
148   /* get the command line options if there are any and set them */
149   PetscCall(TSSetFromOptions(ts));
150 
151 #if defined(PETSC_HAVE_SUNDIALS2)
152   {
153     TSType    type;
154     PetscBool sundialstype = PETSC_FALSE;
155     PetscCall(TSGetType(ts, &type));
156     PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSSUNDIALS, &sundialstype));
157     PetscCheck(!sundialstype || !appctx.useAlhs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use Alhs formulation for TSSUNDIALS type");
158   }
159 #endif
160   /* Sets the initial solution */
161   PetscCall(TSSetSolution(ts, init_sol));
162 
163   stepsz[0] = 1.0 / (2.0 * (nz - 1) * (nz - 1)); /* (mesh_size)^2/2.0 */
164   ftime     = 0.0;
165   for (k = 0; k < nphase; k++) {
166     if (nphase > 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Phase %" PetscInt_FMT " initial time %g, stepsz %g, duration: %g\n", k, (double)ftime, (double)stepsz[k], (double)((k + 1) * T)));
167     PetscCall(TSSetTime(ts, ftime));
168     PetscCall(TSSetTimeStep(ts, stepsz[k]));
169     PetscCall(TSSetMaxSteps(ts, max_steps));
170     PetscCall(TSSetMaxTime(ts, (k + 1) * T));
171     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
172 
173     /* loop over time steps */
174     /*----------------------*/
175     PetscCall(TSSolve(ts, init_sol));
176     PetscCall(TSGetSolveTime(ts, &ftime));
177     PetscCall(TSGetStepNumber(ts, &steps));
178     stepsz[k + 1] = stepsz[k] * 1.5; /* change step size for the next phase */
179   }
180 
181   /* free space */
182   PetscCall(TSDestroy(&ts));
183   PetscCall(MatDestroy(&appctx.Amat));
184   PetscCall(MatDestroy(&Jmat));
185   PetscCall(VecDestroy(&appctx.ksp_rhs));
186   PetscCall(VecDestroy(&appctx.ksp_sol));
187   PetscCall(VecDestroy(&init_sol));
188   PetscCall(VecDestroy(&appctx.solution));
189   PetscCall(PetscFree(z));
190 
191   PetscCall(PetscFinalize());
192   return 0;
193 }
194 
195 /*------------------------------------------------------------------------
196   Set exact solution
197   u(z,t) = sin(6*PI*z)*exp(-36.*PI*PI*t) + 3.*sin(2*PI*z)*exp(-4.*PI*PI*t)
198 --------------------------------------------------------------------------*/
199 PetscScalar exact(PetscScalar z, PetscReal t) {
200   PetscScalar val, ex1, ex2;
201 
202   ex1 = PetscExpReal(-36. * PETSC_PI * PETSC_PI * t);
203   ex2 = PetscExpReal(-4. * PETSC_PI * PETSC_PI * t);
204   val = PetscSinScalar(6 * PETSC_PI * z) * ex1 + 3. * PetscSinScalar(2 * PETSC_PI * z) * ex2;
205   return val;
206 }
207 
208 /*
209    Monitor - User-provided routine to monitor the solution computed at
210    each timestep.  This example plots the solution and computes the
211    error in two different norms.
212 
213    Input Parameters:
214    ts     - the timestep context
215    step   - the count of the current step (with 0 meaning the
216              initial condition)
217    time   - the current time
218    u      - the solution at this timestep
219    ctx    - the user-provided context for this monitoring routine.
220             In this case we use the application context which contains
221             information about the problem size, workspace and the exact
222             solution.
223 */
224 PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) {
225   AppCtx      *appctx = (AppCtx *)ctx;
226   PetscInt     i, m = appctx->m;
227   PetscReal    norm_2, norm_max, h = 1.0 / (m + 1);
228   PetscScalar *u_exact;
229 
230   /* Compute the exact solution */
231   PetscCall(VecGetArrayWrite(appctx->solution, &u_exact));
232   for (i = 0; i < m; i++) u_exact[i] = exact(appctx->z[i + 1], time);
233   PetscCall(VecRestoreArrayWrite(appctx->solution, &u_exact));
234 
235   /* Print debugging information if desired */
236   if (appctx->debug) {
237     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Computed solution vector at time %g\n", (double)time));
238     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
239     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Exact solution vector\n"));
240     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
241   }
242 
243   /* Compute the 2-norm and max-norm of the error */
244   PetscCall(VecAXPY(appctx->solution, -1.0, u));
245   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
246 
247   norm_2 = PetscSqrtReal(h) * norm_2;
248   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
249   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Timestep %" PetscInt_FMT ": time = %g, 2-norm error = %6.4f, max norm error = %6.4f\n", step, (double)time, (double)norm_2, (double)norm_max));
250 
251   /*
252      Print debugging information if desired
253   */
254   if (appctx->debug) {
255     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error vector\n"));
256     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
257   }
258   return 0;
259 }
260 
261 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
262       Function to solve a linear system using KSP
263 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
264 
265 PetscErrorCode Petsc_KSPSolve(AppCtx *obj) {
266   KSP ksp;
267   PC  pc;
268 
269   /*create the ksp context and set the operators,that is, associate the system matrix with it*/
270   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
271   PetscCall(KSPSetOperators(ksp, obj->Amat, obj->Amat));
272 
273   /*get the preconditioner context, set its type and the tolerances*/
274   PetscCall(KSPGetPC(ksp, &pc));
275   PetscCall(PCSetType(pc, PCLU));
276   PetscCall(KSPSetTolerances(ksp, 1.e-7, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
277 
278   /*get the command line options if there are any and set them*/
279   PetscCall(KSPSetFromOptions(ksp));
280 
281   /*get the linear system (ksp) solve*/
282   PetscCall(KSPSolve(ksp, obj->ksp_rhs, obj->ksp_sol));
283 
284   PetscCall(KSPDestroy(&ksp));
285   return 0;
286 }
287 
288 /***********************************************************************
289   Function to return value of basis function or derivative of basis function.
290  ***********************************************************************
291 
292         Arguments:
293           x       = array of xpoints or nodal values
294           xx      = point at which the basis function is to be
295                       evaluated.
296           il      = interval containing xx.
297           iq      = indicates which of the two basis functions in
298                       interval intrvl should be used
299           nll     = array containing the endpoints of each interval.
300           id      = If id ~= 2, the value of the basis function
301                       is calculated; if id = 2, the value of the
302                       derivative of the basis function is returned.
303  ***********************************************************************/
304 
305 PetscScalar bspl(PetscScalar *x, PetscScalar xx, PetscInt il, PetscInt iq, PetscInt nll[][2], PetscInt id) {
306   PetscScalar x1, x2, bfcn;
307   PetscInt    i1, i2, iq1, iq2;
308 
309   /* Determine which basis function in interval intrvl is to be used in */
310   iq1 = iq;
311   if (iq1 == 0) iq2 = 1;
312   else iq2 = 0;
313 
314   /*    Determine endpoint of the interval intrvl   */
315   i1 = nll[il][iq1];
316   i2 = nll[il][iq2];
317 
318   /*   Determine nodal values at the endpoints of the interval intrvl   */
319   x1 = x[i1];
320   x2 = x[i2];
321 
322   /*   Evaluate basis function   */
323   if (id == 2) bfcn = (1.0) / (x1 - x2);
324   else bfcn = (xx - x2) / (x1 - x2);
325   return bfcn;
326 }
327 
328 /*---------------------------------------------------------
329   Function called by rhs function to get B and g
330 ---------------------------------------------------------*/
331 PetscErrorCode femBg(PetscScalar btri[][3], PetscScalar *f, PetscInt nz, PetscScalar *z, PetscReal t) {
332   PetscInt    i, j, jj, il, ip, ipp, ipq, iq, iquad, iqq;
333   PetscInt    nli[num_z][2], indx[num_z];
334   PetscScalar dd, dl, zip, zipq, zz, b_z, bb_z, bij;
335   PetscScalar zquad[num_z][3], dlen[num_z], qdwt[3];
336 
337   /*  initializing everything - btri and f are initialized in rhs.c  */
338   for (i = 0; i < nz; i++) {
339     nli[i][0]   = 0;
340     nli[i][1]   = 0;
341     indx[i]     = 0;
342     zquad[i][0] = 0.0;
343     zquad[i][1] = 0.0;
344     zquad[i][2] = 0.0;
345     dlen[i]     = 0.0;
346   } /*end for (i)*/
347 
348   /*  quadrature weights  */
349   qdwt[0] = 1.0 / 6.0;
350   qdwt[1] = 4.0 / 6.0;
351   qdwt[2] = 1.0 / 6.0;
352 
353   /* 1st and last nodes have Dirichlet boundary condition -
354      set indices there to -1 */
355 
356   for (i = 0; i < nz - 1; i++) indx[i] = i - 1;
357   indx[nz - 1] = -1;
358 
359   ipq = 0;
360   for (il = 0; il < nz - 1; il++) {
361     ip           = ipq;
362     ipq          = ip + 1;
363     zip          = z[ip];
364     zipq         = z[ipq];
365     dl           = zipq - zip;
366     zquad[il][0] = zip;
367     zquad[il][1] = (0.5) * (zip + zipq);
368     zquad[il][2] = zipq;
369     dlen[il]     = PetscAbsScalar(dl);
370     nli[il][0]   = ip;
371     nli[il][1]   = ipq;
372   }
373 
374   for (il = 0; il < nz - 1; il++) {
375     for (iquad = 0; iquad < 3; iquad++) {
376       dd = (dlen[il]) * (qdwt[iquad]);
377       zz = zquad[il][iquad];
378 
379       for (iq = 0; iq < 2; iq++) {
380         ip  = nli[il][iq];
381         b_z = bspl(z, zz, il, iq, nli, 2);
382         i   = indx[ip];
383 
384         if (i > -1) {
385           for (iqq = 0; iqq < 2; iqq++) {
386             ipp  = nli[il][iqq];
387             bb_z = bspl(z, zz, il, iqq, nli, 2);
388             j    = indx[ipp];
389             bij  = -b_z * bb_z;
390 
391             if (j > -1) {
392               jj = 1 + j - i;
393               btri[i][jj] += bij * dd;
394             } else {
395               f[i] += bij * dd * exact(z[ipp], t);
396               /* f[i] += 0.0; */
397               /* if (il==0 && j==-1) { */
398               /* f[i] += bij*dd*exact(zz,t); */
399               /* }*/ /*end if*/
400             }        /*end else*/
401           }          /*end for (iqq)*/
402         }            /*end if (i>0)*/
403       }              /*end for (iq)*/
404     }                /*end for (iquad)*/
405   }                  /*end for (il)*/
406   return 0;
407 }
408 
409 PetscErrorCode femA(AppCtx *obj, PetscInt nz, PetscScalar *z) {
410   PetscInt    i, j, il, ip, ipp, ipq, iq, iquad, iqq;
411   PetscInt    nli[num_z][2], indx[num_z];
412   PetscScalar dd, dl, zip, zipq, zz, bb, bbb, aij;
413   PetscScalar rquad[num_z][3], dlen[num_z], qdwt[3], add_term;
414 
415   /*  initializing everything  */
416   for (i = 0; i < nz; i++) {
417     nli[i][0]   = 0;
418     nli[i][1]   = 0;
419     indx[i]     = 0;
420     rquad[i][0] = 0.0;
421     rquad[i][1] = 0.0;
422     rquad[i][2] = 0.0;
423     dlen[i]     = 0.0;
424   } /*end for (i)*/
425 
426   /*  quadrature weights  */
427   qdwt[0] = 1.0 / 6.0;
428   qdwt[1] = 4.0 / 6.0;
429   qdwt[2] = 1.0 / 6.0;
430 
431   /* 1st and last nodes have Dirichlet boundary condition -
432      set indices there to -1 */
433 
434   for (i = 0; i < nz - 1; i++) indx[i] = i - 1;
435   indx[nz - 1] = -1;
436 
437   ipq = 0;
438 
439   for (il = 0; il < nz - 1; il++) {
440     ip           = ipq;
441     ipq          = ip + 1;
442     zip          = z[ip];
443     zipq         = z[ipq];
444     dl           = zipq - zip;
445     rquad[il][0] = zip;
446     rquad[il][1] = (0.5) * (zip + zipq);
447     rquad[il][2] = zipq;
448     dlen[il]     = PetscAbsScalar(dl);
449     nli[il][0]   = ip;
450     nli[il][1]   = ipq;
451   } /*end for (il)*/
452 
453   for (il = 0; il < nz - 1; il++) {
454     for (iquad = 0; iquad < 3; iquad++) {
455       dd = (dlen[il]) * (qdwt[iquad]);
456       zz = rquad[il][iquad];
457 
458       for (iq = 0; iq < 2; iq++) {
459         ip = nli[il][iq];
460         bb = bspl(z, zz, il, iq, nli, 1);
461         i  = indx[ip];
462         if (i > -1) {
463           for (iqq = 0; iqq < 2; iqq++) {
464             ipp = nli[il][iqq];
465             bbb = bspl(z, zz, il, iqq, nli, 1);
466             j   = indx[ipp];
467             aij = bb * bbb;
468             if (j > -1) {
469               add_term = aij * dd;
470               PetscCall(MatSetValue(obj->Amat, i, j, add_term, ADD_VALUES));
471             } /*endif*/
472           }   /*end for (iqq)*/
473         }     /*end if (i>0)*/
474       }       /*end for (iq)*/
475     }         /*end for (iquad)*/
476   }           /*end for (il)*/
477   PetscCall(MatAssemblyBegin(obj->Amat, MAT_FINAL_ASSEMBLY));
478   PetscCall(MatAssemblyEnd(obj->Amat, MAT_FINAL_ASSEMBLY));
479   return 0;
480 }
481 
482 /*---------------------------------------------------------
483         Function to fill the rhs vector with
484         By + g values ****
485 ---------------------------------------------------------*/
486 PetscErrorCode rhs(AppCtx *obj, PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t) {
487   PetscInt    i, j, js, je, jj;
488   PetscScalar val, g[num_z], btri[num_z][3], add_term;
489 
490   for (i = 0; i < nz - 2; i++) {
491     for (j = 0; j <= 2; j++) btri[i][j] = 0.0;
492     g[i] = 0.0;
493   }
494 
495   /*  call femBg to set the tri-diagonal b matrix and vector g  */
496   femBg(btri, g, nz, z, t);
497 
498   /*  setting the entries of the right hand side vector  */
499   for (i = 0; i < nz - 2; i++) {
500     val = 0.0;
501     js  = 0;
502     if (i == 0) js = 1;
503     je = 2;
504     if (i == nz - 2) je = 1;
505 
506     for (jj = js; jj <= je; jj++) {
507       j = i + jj - 1;
508       val += (btri[i][jj]) * (y[j]);
509     }
510     add_term = val + g[i];
511     PetscCall(VecSetValue(obj->ksp_rhs, (PetscInt)i, (PetscScalar)add_term, INSERT_VALUES));
512   }
513   PetscCall(VecAssemblyBegin(obj->ksp_rhs));
514   PetscCall(VecAssemblyEnd(obj->ksp_rhs));
515   return 0;
516 }
517 
518 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
519 %%   Function to form the right hand side of the time-stepping problem.                       %%
520 %% -------------------------------------------------------------------------------------------%%
521   if (useAlhs):
522     globalout = By+g
523   else if (!useAlhs):
524     globalout = f(y,t)=Ainv(By+g),
525       in which the ksp solver to transform the problem A*ydot=By+g
526       to the problem ydot=f(y,t)=inv(A)*(By+g)
527 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
528 
529 PetscErrorCode RHSfunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) {
530   AppCtx            *obj = (AppCtx *)ctx;
531   PetscScalar        soln[num_z];
532   const PetscScalar *soln_ptr;
533   PetscInt           i, nz = obj->nz;
534   PetscReal          time;
535 
536   /* get the previous solution to compute updated system */
537   PetscCall(VecGetArrayRead(globalin, &soln_ptr));
538   for (i = 0; i < num_z - 2; i++) soln[i] = soln_ptr[i];
539   PetscCall(VecRestoreArrayRead(globalin, &soln_ptr));
540   soln[num_z - 1] = 0.0;
541   soln[num_z - 2] = 0.0;
542 
543   /* clear out the matrix and rhs for ksp to keep things straight */
544   PetscCall(VecSet(obj->ksp_rhs, (PetscScalar)0.0));
545 
546   time = t;
547   /* get the updated system */
548   rhs(obj, soln, nz, obj->z, time); /* setup of the By+g rhs */
549 
550   /* do a ksp solve to get the rhs for the ts problem */
551   if (obj->useAlhs) {
552     /* ksp_sol = ksp_rhs */
553     PetscCall(VecCopy(obj->ksp_rhs, globalout));
554   } else {
555     /* ksp_sol = inv(Amat)*ksp_rhs */
556     PetscCall(Petsc_KSPSolve(obj));
557     PetscCall(VecCopy(obj->ksp_sol, globalout));
558   }
559   return 0;
560 }
561 
562 /*TEST
563 
564     build:
565       requires: !complex
566 
567     test:
568       suffix: euler
569       output_file: output/ex3.out
570 
571     test:
572       suffix: 2
573       args:   -useAlhs
574       output_file: output/ex3.out
575       TODO: Broken because SNESComputeJacobianDefault is incompatible with TSComputeIJacobianConstant
576 
577 TEST*/
578