Actual source code: gbeuler.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: gbeuler.c,v 1.24 2000/01/10 03:18:55 knepley Exp $";
3: #endif
4: /*
5: Code for Timestepping with implicit backwards Euler using Grid Solvers.
6: */
7: #include src/ts/tsimpl.h
8: #include gsolver.h
10: typedef struct {
11: Vec update; /* work vector where new solution is formed */
12: Vec func; /* work vector where F(t[i],u[i]) is stored */
13: Vec rhs; /* work vector for RHS; vec_sol/dt */
14: } TS_BEuler;
16: extern int TSCreate_BEuler(TS);
20: static int GTSStep_BEuler_Linear(GTS ts, int *steps, PetscReal *ltime)
21: {
22: TS_BEuler *beuler = (TS_BEuler*) ts->data;
23: int max_steps = ts->max_steps;
24: Grid grid;
25: MatStructure str;
26: PetscTruth isTimeDependent;
27: PetscScalar mdt;
28: int numFields, field, its;
29: int step, f;
30: int ierr;
33: *steps = -ts->steps;
34: TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
36: /* Set initial guess to be previous solution */
37: VecCopy(ts->vec_sol, beuler->update);
39: /* Step through time */
40: GTSGetGrid(ts, &grid);
41: GridGetNumActiveFields(grid, &numFields);
42: for(step = 0; step < max_steps; step++) {
43: /* Update mesh */
44: ts->time_step_old = ts->time_step;
45: (*ts->ops->update)(ts, ts->ptime, &ts->time_step);
46: if (ts->time_step != ts->time_step_old) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Variable time step not allowed");
48: /* Increment time */
49: ts->ptime += ts->time_step;
50: mdt = 1.0/ts->time_step;
51: if (ts->ptime > ts->max_time) break;
53: /* Update the boundary values */
54: GTSCalcBCValues(ts);
56: /* Compute A(t) and B(t) */
57: GTSEvaluateSystemMatrix(ts, ts->ptime, &ts->A, &ts->B, &str, (PetscObject) ts->jacP);
59: /* Compute f + \Delta t^{-1} I U^n */
60: GTSEvaluateRhs(ts, ts->ptime, ts->vec_sol, beuler->rhs, (PetscObject) ts->funP);
61: for(f = 0; f < numFields; f++) {
62: GridGetActiveField(grid, f, &field);
63: GTSGetTimeDependence(ts, field, &isTimeDependent);
64: if (isTimeDependent == PETSC_TRUE) {
65: GVecEvaluateOperatorGalerkin(beuler->rhs, ts->vec_sol, ts->vec_sol, field, field, IDENTITY, mdt, ts->funP);
66: }
67: }
68: (*ts->ops->applyrhsbc)(ts, beuler->rhs, ts->funP);
70: /* Solve (\Delta t^{-1} I + A) U^{n+1} = \Delta t^{-1} I U^n */
71: SLESSetOperators(ts->sles, ts->A, ts->B, str);
72: SLESSolve(ts->sles, beuler->rhs, beuler->update, &its);
73: ts->linear_its += PetscAbsInt(its);
74: ts->steps++;
76: /* Update solution U^n --> U^{n+1} */
77: VecCopy(beuler->update, ts->vec_sol);
78: TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
79: }
81: *steps += ts->steps;
82: *ltime = ts->ptime;
83: return(0);
84: }
88: static int GTSStep_BEuler_Nonlinear(GTS ts, int *steps, PetscReal *ltime)
89: {
90: TS_BEuler *beuler = (TS_BEuler*) ts->data;
91: Grid grid;
92: int max_steps = ts->max_steps;
93: int its, lits, numFields;
94: int step, field, f;
95: int ierr;
96:
98: *steps = -ts->steps;
99: TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
101: /* Set initial guess to be previous solution */
102: VecCopy(ts->vec_sol, beuler->update);
103: GTSGetGrid(ts, &grid);
104: GridGetNumActiveFields(grid, &numFields);
106: /* Step through time */
107: for(step = 0; step < max_steps; step++)
108: {
109: /* Update mesh */
110: ts->time_step_old = ts->time_step;
111: (*ts->ops->update)(ts, ts->ptime, &ts->time_step);
112: if (ts->time_step != ts->time_step_old)
113: {
114: /* Multiply \Delta t^{-1}_0 I by {\Delta t_0 \over \Delta t} */
115: for(f = 0; f < numFields; f++)
116: {
117: GridGetActiveField(grid, f, &field);
118: if (ts->Iindex[field] >= 0)
119: {GridScaleMatOperator(grid, ts->time_step_old/ts->time_step, ts->Iindex[field]); }
120: }
121: ts->time_step_old = ts->time_step;
122: }
124: /* Increment time */
125: ts->ptime += ts->time_step;
126: if (ts->ptime > ts->max_time) break;
128: /* Update the boundary values */
129: GTSCalcBCValues(ts);
131: /* Solve \Delta t^{-1} I (U^{n+1} - U^n) - F(U^{n+1}, t^{n+1}) = 0 */
132: SNESSolve(ts->snes, beuler->update, &its);
133: SNESGetNumberLinearIterations(ts->snes, &lits);
134: ts->nonlinear_its += PetscAbsInt(its);
135: ts->linear_its += PetscAbsInt(lits);
136: ts->steps++;
138: /* Update solution U^n --> U^{n+1} */
139: VecCopy(beuler->update, ts->vec_sol);
140: TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
141: }
143: *steps += ts->steps;
144: *ltime = ts->ptime;
145: return(0);
146: }
148: /*
149: This defines the nonlinear equation that is to be solved with GSNES
151: \Delta t^{-1} I (U^{n+1} - U^{n}) - F(U^{n+1}, t^{n+1})
152: */
155: int GTSBEulerFunction(GSNES snes, GVec x, GVec y, void *ctx)
156: {
157: Grid grid;
158: PetscObject oldCtx;
159: TS ts;
160: PetscTruth isTimeDependent;
161: PetscScalar mone = -1.0;
162: PetscScalar mdt;
163: int numFields, field, f;
164: int ierr;
167: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
169: mdt = 1.0/ts->time_step;
170: /* Make -F(U^{n+1}, t^{n+1}) */
171: GTSEvaluateRhs(ts, ts->ptime, x, y, (PetscObject) ts->funP);
173: /* Add \Delta t^{-1} I (U^{n+1} - U^n) for each time dependent field */
174: GTSGetGrid(ts, &grid);
175: GTSCreateContext(ts, ts->ptime, (PetscObject) ts->funP, &oldCtx);
176: VecWAXPY(&mone, ts->vec_sol, x, ts->work[0]);
177: GridCalcBCValuesDifference(grid);
178: GridSetBCValuesType(grid, BC_VALUES_DIFF);
179: GridGetNumActiveFields(grid, &numFields);
180: for(f = 0; f < numFields; f++) {
181: GridGetActiveField(grid, f, &field);
182: GTSGetTimeDependence(ts, field, &isTimeDependent);
183: if (isTimeDependent) {
184: GVecEvaluateOperatorGalerkin(y, ts->work[0], ts->work[0], field, field, IDENTITY, mdt, ts->funP);
185: }
186: }
187: GridSetBCValuesType(grid, BC_VALUES);
189: /* Apply boundary conditions */
190: (*ts->ops->applyrhsbc)(ts, y, ts->funP);
192: GTSDestroyContext(ts, (PetscObject) ts->funP, oldCtx);
193: return(0);
194: }
195: /*
196: This defines the nonlinear equation that is to be solved with GSNES
198: P^T (\Delta t^{-1} I P (U^{n+1} - U^{n}) - F(P U^{n+1}, t^{n+1})) + \Delta t^{-1} I_P (U^{n+1}_P - U^n_P) - F_g
199: */
202: int GTSBEulerFunctionConstrained(GSNES snes, GVec x, GVec y, void *ctx)
203: {
204: Grid grid;
205: PetscObject oldCtx;
206: TS ts;
207: Mat P;
208: PetscTruth isTimeDependent;
209: PetscScalar mdt;
210: PetscScalar mone = -1.0;
211: int numFields, field, f;
212: int ierr;
215: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
217: mdt = 1.0/ts->time_step;
218: #ifdef PETSC_USE_BOPT_g
219: PetscTrValid(__LINE__, __FUNCT__, __FILE__, __SDIR__);
220: #endif
222: /* Make -F(P U^{n+1}, t^{n+1}) */
223: GTSGetGrid(ts, &grid);
224: GTSCreateContext(ts, ts->ptime, (PetscObject) ts->funP, &oldCtx);
225: GridGetConstraintMatrix(grid, &P);
226: MatMult(P, x, ts->work[0]);
227: GTSEvaluateRhs(ts, ts->ptime, ts->work[0], ts->work[1], (PetscObject) ts->funP);
229: /* Add \Delta t^{-1} I P (U^{n+1} - U^n) for each time dependent field */
230: MatMult(P, ts->vec_sol, ts->work[2]);
231: VecWAXPY(&mone, ts->work[2], ts->work[0], ts->work[0]);
232: GridCalcBCValuesDifference(grid);
233: GridSetBCValuesType(grid, BC_VALUES_DIFF);
234: GridGetNumActiveFields(grid, &numFields);
235: for(f = 0; f < numFields; f++) {
236: GridGetActiveField(grid, f, &field);
237: GTSGetTimeDependence(ts, field, &isTimeDependent);
238: if (isTimeDependent) {
239: GVecEvaluateOperatorGalerkin(ts->work[1], ts->work[0], ts->work[0], field, field, IDENTITY, mdt, ts->funP);
240: }
241: }
242: GridSetBCValuesType(grid, BC_VALUES);
244: /* Reduce the system with constraints: apply P^T */
245: MatMultTranspose(P, ts->work[1], y);
247: /* Add in extra degrees of freedom */
248: (*((PetscConstraintObject) ctx)->ops->applyrhs)(grid, x, y);
250: /* Apply boundary conditions */
251: (*ts->ops->applyrhsbc)(ts, y, ts->funP);
253: GTSDestroyContext(ts, (PetscObject) ts->funP, oldCtx);
254: #ifdef PETSC_USE_BOPT_g
255: PetscTrValid(__LINE__, __FUNCT__, __FILE__, __SDIR__);
256: #endif
257: return(0);
258: }
260: /*
261: This constructs the Jacobian needed for GSNES
263: J = \Delta t^{-1} I - J_{F}
265: where J_{F} is the given Jacobian of F.
266: */
269: int GTSBEulerJacobian(GSNES snes, GVec x, GMat *J, GMat *M, MatStructure *str, void *ctx)
270: {
271: GTS ts;
275: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
277: GTSEvaluateJacobian(ts, ts->ptime, x, J, M, str, (PetscObject) ctx);
278: return(0);
279: }
281: /*
282: This constructs the Jacobian needed for GSNES
284: J = \Delta t^{-1} I - J_{F}(P U^n)
286: where J_{F} is the given Jacobian of F.
287: */
290: int GTSBEulerJacobianConstrained(GSNES snes, GVec x, GMat *J, GMat *M, MatStructure *str, void *ctx)
291: {
292: GTS ts;
293: Grid grid;
294: Mat P;
295: int ierr;
298: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
300: GTSGetGrid(ts, &grid);
301: GridGetConstraintMatrix(grid, &P);
302: MatMult(P, x, ts->work[0]);
303: GTSEvaluateJacobian(ts, ts->ptime, ts->work[0], J, M, str, (PetscObject) ctx);
304: return(0);
305: }
309: static int GTSSetUp_BEuler_Linear(GTS ts) {
310: TS_BEuler *beuler = (TS_BEuler*) ts->data;
311: Grid grid;
312: int ierr;
315: GTSGetGrid(ts, &grid);
316: VecDuplicate(ts->vec_sol, &beuler->update);
317: VecDuplicate(ts->vec_sol, &beuler->rhs);
318: return(0);
319: }
323: static int GTSSetupGSNES_BEuler(GTS ts)
324: {
325: GVec func;
326: GMat A, B;
327: int ierr;
330: SNESGetFunction(ts->snes, &func, PETSC_NULL, PETSC_NULL);
331: PetscObjectCompose((PetscObject) ts->funP, "TS", (PetscObject) ts);
332: SNESSetFunction(ts->snes, func, GTSBEulerFunction, ts->funP);
333: SNESGetJacobian(ts->snes, &A, &B, PETSC_NULL, PETSC_NULL);
334: PetscObjectCompose((PetscObject) ts->jacP, "TS", (PetscObject) ts);
335: SNESSetJacobian(ts->snes, A, B, GTSBEulerJacobian, ts->jacP);
337: /* This requires that ts is the context for the function */
338: SNESSetSolutionBC(ts->snes, GTSSolutionBCforGSNES);
340: return(0);
341: }
345: static int GTSSetUp_BEuler_Nonlinear(GTS ts)
346: {
347: TS_BEuler *beuler = (TS_BEuler*) ts->data;
348: Grid grid;
349: int ierr;
352: GTSGetGrid(ts, &grid);
354: /* Setup Rhs \Delta t^{-1} I (U^{n+1} - U^{n}) - F(U^{n+1}, t^{n+1}) */
355: VecDuplicate(ts->vec_sol, &beuler->update);
356: ts->nwork = 1;
357: VecDuplicateVecs(ts->vec_sol, ts->nwork, &ts->work);
359: /* Setup the nonlinear solver */
360: GTSSetupGSNES_BEuler(ts);
362: return(0);
363: }
367: static int GTSSetupGSNES_BEuler_Constrained(GTS ts)
368: {
369: GVec func;
370: GMat A, B;
371: int ierr;
374: SNESGetFunction(ts->snes, &func, PETSC_NULL, PETSC_NULL);
375: PetscObjectCompose((PetscObject) ts->funP, "TS", (PetscObject) ts);
376: SNESSetFunction(ts->snes, func, GTSBEulerFunctionConstrained, ts->funP);
377: SNESGetJacobian(ts->snes, &A, &B, PETSC_NULL, PETSC_NULL);
378: PetscObjectCompose((PetscObject) ts->jacP, "TS", (PetscObject) ts);
379: SNESSetJacobian(ts->snes, A, B, GTSBEulerJacobianConstrained, ts->jacP);
381: /* This requires that ts is the context for the function */
382: SNESSetSolutionBC(ts->snes, GTSSolutionBCforGSNES);
384: return(0);
385: }
389: static int GTSSetUp_BEuler_Nonlinear_Constrained(GTS ts)
390: {
391: TS_BEuler *beuler = (TS_BEuler*) ts->data;
392: Grid grid;
393: PetscTruth explicitConst;
394: int ierr;
397: GTSGetGrid(ts, &grid);
399: /* Setup Rhs \Delta t^{-1} I (U^{n+1} - U^{n}) - F(U^{n+1}, t^{n+1}) */
400: VecDuplicate(ts->vec_sol, &beuler->update);
401: ts->nwork = 3;
402: PetscMalloc(ts->nwork * sizeof(GVec *), &ts->work);
404: /* Setup the nonlinear solver */
405: GridGetExplicitConstraints(grid, &explicitConst);
406: if (explicitConst == PETSC_FALSE) {
407: GTSSetupGSNES_BEuler_Constrained(ts);
408: GVecCreate(grid, &ts->work[0]);
409: } else {
410: GTSSetupGSNES_BEuler(ts);
411: GVecCreateConstrained(grid, &ts->work[0]);
412: }
413: GVecDuplicate(ts->work[0], &ts->work[1]);
414: GVecDuplicate(ts->work[0], &ts->work[2]);
416: /* Set the context */
417: GTSCreateConstraintContext(ts);
419: return(0);
420: }
424: static int GTSReform_BEuler(GTS ts)
425: {
427: return(0);
428: }
432: static int GTSReallocate_BEuler(GTS ts)
433: {
434: TS_BEuler *beuler = (TS_BEuler*) ts->data;
435: Grid grid;
436: GSNES newSnes;
437: int ierr;
440: /* Destroy old structures */
441: VecDestroy(beuler->update);
442: if (ts->nwork) {
443: VecDestroyVecs(ts->work, ts->nwork);
444: }
446: /* Recreate GSNES */
447: if (ts->snes) {
448: GTSGetGrid(ts, &grid);
449: GSNESCreate(grid, ts, &newSnes);
450: SNESSetFromOptions(newSnes);
451: GSNESDuplicateMonitors(ts->snes, newSnes);
452: GSNESDestroy(ts->snes);
453: ts->snes = newSnes;
454: }
456: /* Recreate structures */
457: TSSetUp(ts);
458:
459: /* Start off with the current solution */
460: VecCopy(ts->vec_sol, beuler->update);
461: return(0);
462: }
466: static int GTSDestroy_BEuler(TS ts)
467: {
468: TS_BEuler *beuler = (TS_BEuler*) ts->data;
469: int ierr;
472: if (beuler->update) {
473: VecDestroy(beuler->update);
474: }
475: if (beuler->func) {
476: VecDestroy(beuler->func);
477: }
478: if (beuler->rhs) {
479: VecDestroy(beuler->rhs);
480: }
481: if (ts->nwork) {
482: VecDestroyVecs(ts->work, ts->nwork);
483: }
484: PetscFree(beuler);
485: return(0);
486: }
490: static int GTSSetFromOptions_BEuler_Linear(TS ts)
491: {
493: return(0);
494: }
498: static int GTSSetFromOptions_BEuler_Nonlinear(TS ts)
499: {
501: return(0);
502: }
506: static int GTSPrintHelp_BEuler(TS ts, char *prefix)
507: {
509: return(0);
510: }
514: static int GTSView_BEuler(TS ts, PetscViewer viewer)
515: {
517: return(0);
518: }
520: EXTERN_C_BEGIN
523: int GTSCreate_BEuler(TS ts)
524: {
525: TS_BEuler *beuler;
526: KSP ksp;
527: Grid grid;
528: PetscTruth flag;
529: int ierr;
532: ts->ops->destroy = GTSDestroy_BEuler;
533: ts->ops->view = GTSView_BEuler;
534: ts->ops->printhelp = GTSPrintHelp_BEuler;
535: PetscObjectChangeSerializeName((PetscObject) ts, GTS_SER_BEULER_BINARY);
537: if (ts->problem_type == TS_LINEAR) {
538: if (!ts->A) SETERRQ(PETSC_ERR_ARG_WRONG, "Must set rhs matrix for linear problem");
539: ts->ops->setup = GTSSetUp_BEuler_Linear;
540: ts->ops->step = GTSStep_BEuler_Linear;
541: ts->ops->setfromoptions = GTSSetFromOptions_BEuler_Linear;
542: ts->ops->reallocate = GTSReallocate_BEuler;
543: ts->ops->reform = GTSReform_BEuler;
544: SLESCreate(ts->comm, &ts->sles);
545: SLESGetKSP(ts->sles, &ksp);
546: KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
547: GTSGetGrid(ts, &grid);
548: GridIsConstrained(grid, &flag);
549: if (flag == PETSC_TRUE) SETERRQ(PETSC_ERR_SUP, "Constraints not supported for the linear problem yet");
550: } else if (ts->problem_type == TS_NONLINEAR) {
551: GTSGetGrid(ts, &grid);
552: GridIsConstrained(grid, &flag);
553: GSNESCreate(grid, ts, &ts->snes);
554: if (flag) {
555: ts->ops->setup = GTSSetUp_BEuler_Nonlinear_Constrained;
556: ts->ops->step = GTSStep_BEuler_Nonlinear;
557: ts->ops->setfromoptions = GTSSetFromOptions_BEuler_Nonlinear;
558: ts->ops->reform = GTSReform_BEuler;
559: ts->ops->reallocate = GTSReallocate_BEuler;
560: } else {
561: ts->ops->setup = GTSSetUp_BEuler_Nonlinear;
562: ts->ops->step = GTSStep_BEuler_Nonlinear;
563: ts->ops->setfromoptions = GTSSetFromOptions_BEuler_Nonlinear;
564: ts->ops->reform = GTSReform_BEuler;
565: ts->ops->reallocate = GTSReallocate_BEuler;
566: }
567: } else {
568: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid problem type");
569: }
571: PetscNew(TS_BEuler, &beuler);
572: PetscLogObjectMemory(ts, sizeof(TS_BEuler));
573: PetscMemzero(beuler, sizeof(TS_BEuler));
574: ts->data = (void *) beuler;
576: return(0);
577: }
578: EXTERN_C_END
580: EXTERN_C_BEGIN
583: int GTSSerialize_BEuler(MPI_Comm comm, TS *ts, PetscViewer viewer, PetscTruth store)
584: {
585: TS t;
586: Grid grid;
587: int fd;
588: int numFields;
589: int hasPC;
590: int zero = 0;
591: int one = 0;
592: int ierr;
595: PetscViewerBinaryGetDescriptor(viewer, &fd);
596: if (store) {
597: t = *ts;
598: PetscBinaryWrite(fd, &t->problem_type, 1, PETSC_INT, 0);
599: PetscBinaryWrite(fd, &t->isGTS, 1, PETSC_INT, 0);
600: GTSGetGrid(t, &grid);
601: GridGetNumFields(grid, &numFields);
602: PetscBinaryWrite(fd, t->isExplicit, numFields, PETSC_INT, 0);
603: PetscBinaryWrite(fd, t->Iindex, numFields, PETSC_INT, 0);
604: if (t->problem_type == TS_LINEAR) {
605: GMatSerialize(grid, &t->A, viewer, store);
606: if (t->B == t->A) {
607: PetscBinaryWrite(fd, &zero, 1, PETSC_INT, 0);
608: } else {
609: PetscBinaryWrite(fd, &one, 1, PETSC_INT, 0);
610: GMatSerialize(grid, &t->B, viewer, store);
611: }
612: }
613: PetscBinaryWrite(fd, &t->max_steps, 1, PETSC_INT, 0);
614: PetscBinaryWrite(fd, &t->max_time, 1, PETSC_SCALAR, 0);
615: PetscBinaryWrite(fd, &t->time_step, 1, PETSC_SCALAR, 0);
616: PetscBinaryWrite(fd, &t->initial_time_step, 1, PETSC_SCALAR, 0);
617: PetscBinaryWrite(fd, &t->max_time, 1, PETSC_SCALAR, 0);
618: PetscBinaryWrite(fd, &t->steps, 1, PETSC_INT, 0);
619: PetscBinaryWrite(fd, &t->ptime, 1, PETSC_SCALAR, 0);
620: PetscBinaryWrite(fd, &t->linear_its, 1, PETSC_INT, 0);
621: PetscBinaryWrite(fd, &t->nonlinear_its, 1, PETSC_INT, 0);
623: /* Protect against nasty option overwrites */
624: PetscOptionsClearValue("-ts_init_time");
625: PetscOptionsClearValue("-ts_init_time_step");
626: } else {
627: TSCreate(comm, &t);
629: PetscBinaryRead(fd, &t->problem_type, 1, PETSC_INT);
630: PetscBinaryRead(fd, &t->isGTS, 1, PETSC_INT);
631: t->bops->destroy = (int (*)(PetscObject)) GTSDestroy;
632: t->bops->view = (int (*)(PetscObject, PetscViewer)) GTSView;
633: PetscStrallocpy("gbeuler", &t->type_name);
635: grid = (Grid) *ts;
636: GridGetNumFields(grid, &numFields);
637: PetscMalloc(numFields * sizeof(PetscTruth), &t->isExplicit);
638: PetscBinaryRead(fd, t->isExplicit, numFields, PETSC_INT);
639: PetscMalloc(numFields * sizeof(int), &t->Iindex);
640: PetscBinaryRead(fd, t->Iindex, numFields, PETSC_INT);
641: if (t->problem_type == TS_LINEAR) {
642: GMatSerialize(grid, &t->A, viewer, store);
643: PetscBinaryRead(fd, &hasPC, 1, PETSC_INT);
644: if (hasPC) {
645: GMatSerialize(grid, &t->B, viewer, store);
646: } else {
647: t->B = t->A;
648: }
649: }
650: /* We must have a Grid parent for the constructor */
651: PetscObjectCompose((PetscObject) t, "Grid", (PetscObject) grid);
652: PetscBinaryRead(fd, &t->max_steps, 1, PETSC_INT);
653: PetscBinaryRead(fd, &t->max_time, 1, PETSC_SCALAR);
654: PetscBinaryRead(fd, &t->time_step, 1, PETSC_SCALAR);
655: PetscBinaryRead(fd, &t->initial_time_step, 1, PETSC_SCALAR);
656: PetscBinaryRead(fd, &t->max_time, 1, PETSC_SCALAR);
657: PetscBinaryRead(fd, &t->steps, 1, PETSC_INT);
658: PetscBinaryRead(fd, &t->ptime, 1, PETSC_SCALAR);
659: PetscBinaryRead(fd, &t->linear_its, 1, PETSC_INT);
660: PetscBinaryRead(fd, &t->nonlinear_its, 1, PETSC_INT);
662: GTSCreate_BEuler(t);
663: *ts = t;
664: }
666: return(0);
667: }
668: EXTERN_C_END