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);
18: static int GTSStep_BEuler_Linear(GTS ts, int *steps, PetscReal *ltime)
19: {
20: TS_BEuler *beuler = (TS_BEuler*) ts->data;
21: int max_steps = ts->max_steps;
22: MatStructure str;
23: int its;
24: int step;
25: int ierr;
28: *steps = -ts->steps;
29: ierr = TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
31: /* Set initial guess to be previous solution */
32: VecCopy(ts->vec_sol, beuler->update);
34: /* Step through time */
35: for(step = 0; step < max_steps; step++) {
36: /* Update mesh */
37: ts->time_step_old = ts->time_step;
38: (*ts->ops->update)(ts, ts->ptime, &ts->time_step);
39: if (ts->time_step != ts->time_step_old) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Variable time step not allowed");
41: /* Increment time */
42: ts->ptime += ts->time_step;
43: if (ts->ptime > ts->max_time) break;
45: /* Update the boundary values */
46: GTSCalcBCValues(ts);
48: /* Compute A(t) and B(t) */
49: GTSEvaluateSystemMatrix(ts, ts->ptime, &ts->A, &ts->B, &str, (PetscObject) ts->jacP);
51: /* Compute Delta t^{-1} I U^n */
52: GTSEvaluateRhs(ts, ts->ptime, ts->vec_sol, beuler->rhs, (PetscObject) ts->funP);
53: (*ts->ops->applyrhsbc)(ts, beuler->rhs, ts->funP);
55: /* Solve (Delta t^{-1} I - A) U^{n+1} = Delta t^{-1} I U^n */
56: SLESSetOperators(ts->sles, ts->A, ts->B, str);
57: SLESSolve(ts->sles, beuler->rhs, beuler->update, &its);
58: ts->linear_its += PetscAbsInt(its);
59: ts->steps++;
61: /* Update solution U^n --> U^{n+1} */
62: VecCopy(beuler->update, ts->vec_sol);
63: TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
64: }
66: *steps += ts->steps;
67: *ltime = ts->ptime;
68: return(0);
69: }
71: static int GTSStep_BEuler_Nonlinear(GTS ts, int *steps, PetscReal *ltime)
72: {
73: TS_BEuler *beuler = (TS_BEuler*) ts->data;
74: Grid grid;
75: int max_steps = ts->max_steps;
76: int its, lits, numFields;
77: int step, field, f;
78: int ierr;
79:
81: *steps = -ts->steps;
82: ierr = TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
84: /* Set initial guess to be previous solution */
85: VecCopy(ts->vec_sol, beuler->update);
86: GTSGetGrid(ts, &grid);
87: GridGetNumActiveFields(grid, &numFields);
89: /* Step through time */
90: for(step = 0; step < max_steps; step++)
91: {
92: /* Update mesh */
93: ts->time_step_old = ts->time_step;
94: (*ts->ops->update)(ts, ts->ptime, &ts->time_step);
95: if (ts->time_step != ts->time_step_old)
96: {
97: /* Multiply Delta t^{-1}_0 I by {Delta t_0 over Delta t} */
98: for(f = 0; f < numFields; f++)
99: {
100: GridGetActiveField(grid, f, &field);
101: if (ts->Iindex[field] >= 0)
102: {GridScaleMatOperator(grid, ts->time_step_old/ts->time_step, ts->Iindex[field]); }
103: }
104: ts->time_step_old = ts->time_step;
105: }
107: /* Increment time */
108: ts->ptime += ts->time_step;
109: if (ts->ptime > ts->max_time) break;
111: /* Update the boundary values */
112: GTSCalcBCValues(ts);
114: /* Solve Delta t^{-1} I (U^{n+1} - U^n) - F(U^{n+1}, t^{n+1}) = 0 */
115: SNESSolve(ts->snes, beuler->update, &its);
116: SNESGetNumberLinearIterations(ts->snes, &lits);
117: ts->nonlinear_its += PetscAbsInt(its);
118: ts->linear_its += PetscAbsInt(lits);
119: ts->steps++;
121: /* Update solution U^n --> U^{n+1} */
122: VecCopy(beuler->update, ts->vec_sol);
123: TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol);
124: }
126: *steps += ts->steps;
127: *ltime = ts->ptime;
128: return(0);
129: }
131: /*
132: This defines the nonlinear equation that is to be solved with GSNES
134: Delta t^{-1} I (U^{n+1} - U^{n}) - F(U^{n+1}, t^{n+1})
135: */
136: int GTSBEulerFunction(GSNES snes, GVec x, GVec y, void *ctx)
137: {
138: Grid grid;
139: PetscObject oldCtx;
140: TS ts;
141: PetscTruth isTimeDependent;
142: PetscScalar mone = -1.0;
143: PetscScalar mdt;
144: int numFields, field, f;
145: int ierr;
148: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
150: mdt = 1.0/ts->time_step;
151: /* Make -F(U^{n+1}, t^{n+1}) */
152: GTSEvaluateRhs(ts, ts->ptime, x, y, (PetscObject) ts->funP);
154: /* Add Delta t^{-1} I (U^{n+1} - U^n) for each time dependent field */
155: GTSGetGrid(ts, &grid);
156: GTSCreateContext(ts, ts->ptime, (PetscObject) ts->funP, &oldCtx);
157: VecWAXPY(&mone, ts->vec_sol, x, ts->work[0]);
158: GridCalcBCValuesDifference(grid);
159: GridSetBCValuesType(grid, BC_VALUES_DIFF);
160: GridGetNumActiveFields(grid, &numFields);
161: for(f = 0; f < numFields; f++) {
162: GridGetActiveField(grid, f, &field);
163: GTSGetTimeDependence(ts, field, &isTimeDependent);
164: if (isTimeDependent) {
165: GVecEvaluateOperatorGalerkin(y, ts->work[0], ts->work[0], field, field, IDENTITY, mdt, ts->funP);
166: }
167: }
168: GridSetBCValuesType(grid, BC_VALUES);
170: /* Apply boundary conditions */
171: (*ts->ops->applyrhsbc)(ts, y, ts->funP);
173: GTSDestroyContext(ts, (PetscObject) ts->funP, oldCtx);
174: return(0);
175: }
176: /*
177: This defines the nonlinear equation that is to be solved with GSNES
179: 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
180: */
181: int GTSBEulerFunctionConstrained(GSNES snes, GVec x, GVec y, void *ctx)
182: {
183: Grid grid;
184: PetscObject oldCtx;
185: TS ts;
186: Mat P;
187: PetscTruth isTimeDependent;
188: PetscScalar mdt;
189: PetscScalar mone = -1.0;
190: int numFields, field, f;
191: int ierr;
194: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
196: mdt = 1.0/ts->time_step;
197: #ifdef PETSC_USE_BOPT_g
198: #endif
200: /* Make -F(P U^{n+1}, t^{n+1}) */
201: GTSGetGrid(ts, &grid);
202: GTSCreateContext(ts, ts->ptime, (PetscObject) ts->funP, &oldCtx);
203: GridGetConstraintMatrix(grid, &P);
204: MatMult(P, x, ts->work[0]);
205: GTSEvaluateRhs(ts, ts->ptime, ts->work[0], ts->work[1], (PetscObject) ts->funP);
207: /* Add Delta t^{-1} I P (U^{n+1} - U^n) for each time dependent field */
208: MatMult(P, ts->vec_sol, ts->work[2]);
209: VecWAXPY(&mone, ts->work[2], ts->work[0], ts->work[0]);
210: GridCalcBCValuesDifference(grid);
211: GridSetBCValuesType(grid, BC_VALUES_DIFF);
212: GridGetNumActiveFields(grid, &numFields);
213: for(f = 0; f < numFields; f++) {
214: GridGetActiveField(grid, f, &field);
215: GTSGetTimeDependence(ts, field, &isTimeDependent);
216: if (isTimeDependent) {
217: GVecEvaluateOperatorGalerkin(ts->work[1], ts->work[0], ts->work[0], field, field, IDENTITY, mdt, ts->funP);
218: }
219: }
220: GridSetBCValuesType(grid, BC_VALUES);
222: /* Reduce the system with constraints: apply P^T */
223: MatMultTranspose(P, ts->work[1], y);
225: /* Add in extra degrees of freedom */
226: (*((PetscConstraintObject) ctx)->ops->applyrhs)(grid, x, y);
228: /* Apply boundary conditions */
229: (*ts->ops->applyrhsbc)(ts, y, ts->funP);
231: GTSDestroyContext(ts, (PetscObject) ts->funP, oldCtx);
232: #ifdef PETSC_USE_BOPT_g
233: #endif
234: return(0);
235: }
237: /*
238: This constructs the Jacobian needed for GSNES
240: J = Delta t^{-1} I - J_{F}
242: where J_{F} is the given Jacobian of F.
243: */
244: int GTSBEulerJacobian(GSNES snes, GVec x, GMat *J, GMat *M, MatStructure *str, void *ctx)
245: {
246: GTS ts;
250: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
252: GTSEvaluateJacobian(ts, ts->ptime, x, J, M, str, (PetscObject) ctx);
253: return(0);
254: }
256: /*
257: This constructs the Jacobian needed for GSNES
259: J = Delta t^{-1} I - J_{F}(P U^n)
261: where J_{F} is the given Jacobian of F.
262: */
263: int GTSBEulerJacobianConstrained(GSNES snes, GVec x, GMat *J, GMat *M, MatStructure *str, void *ctx)
264: {
265: GTS ts;
266: Grid grid;
267: Mat P;
268: int ierr;
271: PetscObjectQuery((PetscObject) ctx, "TS", (PetscObject *) &ts);
273: GTSGetGrid(ts, &grid);
274: GridGetConstraintMatrix(grid, &P);
275: MatMult(P, x, ts->work[0]);
276: GTSEvaluateJacobian(ts, ts->ptime, ts->work[0], J, M, str, (PetscObject) ctx);
277: return(0);
278: }
280: #if 0
281: static int GTSSetUp_BEuler_Linear_Constant_Matrix(TS ts)
282: {
283: PetscScalar mdt = 1.0/ts->time_step;
284: Grid grid;
285: PetscTruth isTimeDependent;
286: MatStructure str;
287: int numFields, field, f;
288: int ierr;
291: GTSGetGrid(ts, &grid);
293: /* Setup Rhs Delta t^{-1} I U^{n} */
294: /* Setup system matrix A = Delta t^{-1} I - A_{orig} */
295: GridScaleSystemMatrix(grid, -1.0);
296: GridGetNumFields(grid, &numFields);
297: GridGetNumActiveFields(grid, &numFields);
298: for(f = 0; f < numFields; f++) {
299: GridGetActiveField(grid, f, &field);
300: GTSGetTimeDependence(ts, field, &isTimeDependent);
301: if (isTimeDependent) {
302: GridAddRhsOperator(grid, field, field, IDENTITY, mdt, PETSC_FALSE, PETSC_NULL);
303: GridAddMatOperator(grid, field, field, IDENTITY, mdt, PETSC_FALSE, PETSC_NULL);
304: }
305: }
307: /* Build system matrix and apply boundary conditions */
308: (*ts->ops->rhsmatrix)(ts, ts->ptime, &ts->A, &ts->B, &str, ts->jacP);
309: (*ts->ops->applymatrixbc)(ts, ts->A, ts->B, ts->jacP);
310: SLESSetOperators(ts->sles, ts->A, ts->B, SAME_NONZERO_PATTERN);
312: return(0);
313: }
314: #endif
316: static int GTSSetUp_BEuler_Linear(GTS ts)
317: {
318: TS_BEuler *beuler = (TS_BEuler*) ts->data;
319: PetscScalar mdt = 1.0/ts->time_step;
320: Grid grid;
321: PetscTruth isTimeDependent;
322: int numFields, field, f;
323: int ierr;
326: GTSGetGrid(ts, &grid);
327: VecDuplicate(ts->vec_sol, &beuler->update);
328: VecDuplicate(ts->vec_sol, &beuler->rhs);
330: /* Setup Rhs Delta t^{-1} I U^{n} */
331: /* Setup system matrix A = Delta t^{-1} I - A_{orig} */
332: GridScaleSystemMatrix(grid, -1.0);
333: GridGetNumActiveFields(grid, &numFields);
334: for(f = 0; f < numFields; f++) {
335: GridGetActiveField(grid, f, &field);
336: GTSGetTimeDependence(ts, field, &isTimeDependent);
337: if (isTimeDependent) {
338: GridAddRhsOperator(grid, field, field, IDENTITY, mdt, PETSC_FALSE, PETSC_NULL);
339: GridAddMatOperator(grid, field, field, IDENTITY, mdt, PETSC_FALSE, PETSC_NULL);
340: }
341: }
343: return(0);
344: }
346: static int GTSSetupGSNES_BEuler(GTS ts)
347: {
348: GVec func;
349: GMat A, B;
350: int ierr;
353: SNESGetFunction(ts->snes, &func, PETSC_NULL, PETSC_NULL);
354: PetscObjectCompose((PetscObject) ts->funP, "TS", (PetscObject) ts);
355: SNESSetFunction(ts->snes, func, GTSBEulerFunction, ts->funP);
356: SNESGetJacobian(ts->snes, &A, &B, PETSC_NULL, PETSC_NULL);
357: PetscObjectCompose((PetscObject) ts->jacP, "TS", (PetscObject) ts);
358: SNESSetJacobian(ts->snes, A, B, GTSBEulerJacobian, ts->jacP);
360: /* This requires that ts is the context for the function */
361: SNESSetSolutionBC(ts->snes, GTSSolutionBCforGSNES);
363: return(0);
364: }
366: static int GTSSetUp_BEuler_Nonlinear(GTS ts)
367: {
368: TS_BEuler *beuler = (TS_BEuler*) ts->data;
369: Grid grid;
370: int ierr;
373: GTSGetGrid(ts, &grid);
375: /* Setup Rhs Delta t^{-1} I (U^{n+1} - U^{n}) - F(U^{n+1}, t^{n+1}) */
376: VecDuplicate(ts->vec_sol, &beuler->update);
377: ts->nwork = 1;
378: VecDuplicateVecs(ts->vec_sol, ts->nwork, &ts->work);
380: /* Setup the nonlinear solver */
381: GTSSetupGSNES_BEuler(ts);
383: return(0);
384: }
386: static int GTSSetupGSNES_BEuler_Constrained(GTS ts)
387: {
388: GVec func;
389: GMat A, B;
390: int ierr;
393: SNESGetFunction(ts->snes, &func, PETSC_NULL, PETSC_NULL);
394: PetscObjectCompose((PetscObject) ts->funP, "TS", (PetscObject) ts);
395: SNESSetFunction(ts->snes, func, GTSBEulerFunctionConstrained, ts->funP);
396: SNESGetJacobian(ts->snes, &A, &B, PETSC_NULL, PETSC_NULL);
397: PetscObjectCompose((PetscObject) ts->jacP, "TS", (PetscObject) ts);
398: SNESSetJacobian(ts->snes, A, B, GTSBEulerJacobianConstrained, ts->jacP);
400: /* This requires that ts is the context for the function */
401: SNESSetSolutionBC(ts->snes, GTSSolutionBCforGSNES);
403: return(0);
404: }
406: static int GTSSetUp_BEuler_Nonlinear_Constrained(GTS ts)
407: {
408: TS_BEuler *beuler = (TS_BEuler*) ts->data;
409: Grid grid;
410: PetscTruth explicitConst;
411: int ierr;
414: GTSGetGrid(ts, &grid);
416: /* Setup Rhs Delta t^{-1} I (U^{n+1} - U^{n}) - F(U^{n+1}, t^{n+1}) */
417: VecDuplicate(ts->vec_sol, &beuler->update);
418: ts->nwork = 3;
419: PetscMalloc(ts->nwork * sizeof(GVec *), &ts->work);
421: /* Setup the nonlinear solver */
422: GridGetExplicitConstraints(grid, &explicitConst);
423: if (explicitConst == PETSC_FALSE) {
424: GTSSetupGSNES_BEuler_Constrained(ts);
425: GVecCreate(grid, &ts->work[0]);
426: } else {
427: GTSSetupGSNES_BEuler(ts);
428: GVecCreateConstrained(grid, &ts->work[0]);
429: }
430: GVecDuplicate(ts->work[0], &ts->work[1]);
431: GVecDuplicate(ts->work[0], &ts->work[2]);
433: /* Set the context */
434: GTSCreateConstraintContext(ts);
436: return(0);
437: }
439: static int GTSReform_BEuler(GTS ts)
440: {
442: return(0);
443: }
445: static int GTSReallocate_BEuler(GTS ts)
446: {
447: TS_BEuler *beuler = (TS_BEuler*) ts->data;
448: Grid grid;
449: GSNES newSnes;
450: int ierr;
453: /* Destroy old structures */
454: VecDestroy(beuler->update);
455: if (ts->nwork) {
456: VecDestroyVecs(ts->work, ts->nwork);
457: }
459: /* Recreate GSNES */
460: if (ts->snes) {
461: GTSGetGrid(ts, &grid);
462: GSNESCreate(grid, ts, &newSnes);
463: SNESSetFromOptions(newSnes);
464: GSNESDuplicateMonitors(ts->snes, newSnes);
465: GSNESDestroy(ts->snes);
466: ts->snes = newSnes;
467: }
469: /* Recreate structures */
470: TSSetUp(ts);
471:
472: /* Start off with the current solution */
473: VecCopy(ts->vec_sol, beuler->update);
474: return(0);
475: }
477: static int GTSDestroy_BEuler(TS ts)
478: {
479: TS_BEuler *beuler = (TS_BEuler*) ts->data;
480: int ierr;
483: if (beuler->update) {
484: VecDestroy(beuler->update);
485: }
486: if (beuler->func) {
487: VecDestroy(beuler->func);
488: }
489: if (beuler->rhs) {
490: VecDestroy(beuler->rhs);
491: }
492: if (ts->nwork) {
493: VecDestroyVecs(ts->work, ts->nwork);
494: }
495: PetscFree(beuler);
496: return(0);
497: }
499: static int GTSSetFromOptions_BEuler_Linear(TS ts)
500: {
502: return(0);
503: }
505: static int GTSSetFromOptions_BEuler_Nonlinear(TS ts)
506: {
508: return(0);
509: }
511: static int GTSPrintHelp_BEuler(TS ts, char *prefix)
512: {
514: return(0);
515: }
517: static int GTSView_BEuler(TS ts, PetscViewer viewer)
518: {
520: return(0);
521: }
523: EXTERN_C_BEGIN
524: int GTSCreate_BEuler(TS ts)
525: {
526: TS_BEuler *beuler;
527: KSP ksp;
528: Grid grid;
529: PetscTruth flag;
530: int ierr;
533: ts->ops->destroy = GTSDestroy_BEuler;
534: ts->ops->view = GTSView_BEuler;
535: ts->ops->printhelp = GTSPrintHelp_BEuler;
536: PetscObjectChangeSerializeName((PetscObject) ts, GTS_SER_BEULER_BINARY);
538: if (ts->problem_type == TS_LINEAR) {
539: if (!ts->A) SETERRQ(PETSC_ERR_ARG_WRONG, "Must set rhs matrix for linear problem");
540: ts->ops->setup = GTSSetUp_BEuler_Linear;
541: ts->ops->step = GTSStep_BEuler_Linear;
542: ts->ops->setfromoptions = GTSSetFromOptions_BEuler_Linear;
543: ts->ops->reallocate = GTSReallocate_BEuler;
544: ts->ops->reform = GTSReform_BEuler;
545: SLESCreate(ts->comm, &ts->sles);
546: SLESGetKSP(ts->sles, &ksp);
547: KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
548: GTSGetGrid(ts, &grid);
549: GridIsConstrained(grid, &flag);
550: if (flag == PETSC_TRUE) SETERRQ(PETSC_ERR_SUP, "Constraints not supported for the linear problem yet");
551: } else if (ts->problem_type == TS_NONLINEAR) {
552: GTSGetGrid(ts, &grid);
553: GridIsConstrained(grid, &flag);
554: GSNESCreate(grid, ts, &ts->snes);
555: if (flag) {
556: ts->ops->setup = GTSSetUp_BEuler_Nonlinear_Constrained;
557: ts->ops->step = GTSStep_BEuler_Nonlinear;
558: ts->ops->setfromoptions = GTSSetFromOptions_BEuler_Nonlinear;
559: ts->ops->reform = GTSReform_BEuler;
560: ts->ops->reallocate = GTSReallocate_BEuler;
561: } else {
562: ts->ops->setup = GTSSetUp_BEuler_Nonlinear;
563: ts->ops->step = GTSStep_BEuler_Nonlinear;
564: ts->ops->setfromoptions = GTSSetFromOptions_BEuler_Nonlinear;
565: ts->ops->reform = GTSReform_BEuler;
566: ts->ops->reallocate = GTSReallocate_BEuler;
567: }
568: } else {
569: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid problem type");
570: }
572: PetscNew(TS_BEuler, &beuler);
573: PetscLogObjectMemory(ts, sizeof(TS_BEuler));
574: PetscMemzero(beuler, sizeof(TS_BEuler));
575: ts->data = (void *) beuler;
577: return(0);
578: }
579: EXTERN_C_END
581: EXTERN_C_BEGIN
582: int GTSSerialize_BEuler(MPI_Comm comm, TS *ts, PetscViewer viewer, PetscTruth store)
583: {
584: TS t;
585: Grid grid;
586: int fd;
587: int numFields;
588: int hasPC;
589: int zero = 0;
590: int one = 0;
591: int ierr;
594: PetscViewerBinaryGetDescriptor(viewer, &fd);
595: if (store) {
596: t = *ts;
597: PetscBinaryWrite(fd, &t->problem_type, 1, PETSC_INT, 0);
598: PetscBinaryWrite(fd, &t->isGTS, 1, PETSC_INT, 0);
599: GTSGetGrid(t, &grid);
600: GridGetNumFields(grid, &numFields);
601: PetscBinaryWrite(fd, t->isExplicit, numFields, PETSC_INT, 0);
602: PetscBinaryWrite(fd, t->Iindex, numFields, PETSC_INT, 0);
603: if (t->problem_type == TS_LINEAR) {
604: GMatSerialize(grid, &t->A, viewer, store);
605: if (t->B == t->A) {
606: PetscBinaryWrite(fd, &zero, 1, PETSC_INT, 0);
607: } else {
608: PetscBinaryWrite(fd, &one, 1, PETSC_INT, 0);
609: GMatSerialize(grid, &t->B, viewer, store);
610: }
611: }
612: PetscBinaryWrite(fd, &t->max_steps, 1, PETSC_INT, 0);
613: PetscBinaryWrite(fd, &t->max_time, 1, PETSC_SCALAR, 0);
614: PetscBinaryWrite(fd, &t->time_step, 1, PETSC_SCALAR, 0);
615: PetscBinaryWrite(fd, &t->initial_time_step, 1, PETSC_SCALAR, 0);
616: PetscBinaryWrite(fd, &t->max_time, 1, PETSC_SCALAR, 0);
617: PetscBinaryWrite(fd, &t->steps, 1, PETSC_INT, 0);
618: PetscBinaryWrite(fd, &t->ptime, 1, PETSC_SCALAR, 0);
619: PetscBinaryWrite(fd, &t->linear_its, 1, PETSC_INT, 0);
620: PetscBinaryWrite(fd, &t->nonlinear_its, 1, PETSC_INT, 0);
622: /* Protect against nasty option overwrites */
623: PetscOptionsClearValue("-ts_init_time");
624: PetscOptionsClearValue("-ts_init_time_step");
625: } else {
626: TSCreate(comm, &t);
628: PetscBinaryRead(fd, &t->problem_type, 1, PETSC_INT);
629: PetscBinaryRead(fd, &t->isGTS, 1, PETSC_INT);
630: t->bops->destroy = (int (*)(PetscObject)) GTSDestroy;
631: t->bops->view = (int (*)(PetscObject, PetscViewer)) GTSView;
632: PetscStrallocpy("gbeuler", &t->type_name);
634: grid = (Grid) *ts;
635: GridGetNumFields(grid, &numFields);
636: PetscMalloc(numFields * sizeof(PetscTruth), &t->isExplicit);
637: PetscBinaryRead(fd, t->isExplicit, numFields, PETSC_INT);
638: PetscMalloc(numFields * sizeof(int), &t->Iindex);
639: PetscBinaryRead(fd, t->Iindex, numFields, PETSC_INT);
640: if (t->problem_type == TS_LINEAR) {
641: GMatSerialize(grid, &t->A, viewer, store);
642: PetscBinaryRead(fd, &hasPC, 1, PETSC_INT);
643: if (hasPC) {
644: GMatSerialize(grid, &t->B, viewer, store);
645: } else {
646: t->B = t->A;
647: }
648: }
649: /* We must have a Grid parent for the constructor */
650: PetscObjectCompose((PetscObject) t, "Grid", (PetscObject) grid);
651: PetscBinaryRead(fd, &t->max_steps, 1, PETSC_INT);
652: PetscBinaryRead(fd, &t->max_time, 1, PETSC_SCALAR);
653: PetscBinaryRead(fd, &t->time_step, 1, PETSC_SCALAR);
654: PetscBinaryRead(fd, &t->initial_time_step, 1, PETSC_SCALAR);
655: PetscBinaryRead(fd, &t->max_time, 1, PETSC_SCALAR);
656: PetscBinaryRead(fd, &t->steps, 1, PETSC_INT);
657: PetscBinaryRead(fd, &t->ptime, 1, PETSC_SCALAR);
658: PetscBinaryRead(fd, &t->linear_its, 1, PETSC_INT);
659: PetscBinaryRead(fd, &t->nonlinear_its, 1, PETSC_INT);
661: GTSCreate_BEuler(t);
662: *ts = t;
663: }
665: return(0);
666: }
667: EXTERN_C_END