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