1dd8e379bSPierre Jolivet static char help[] = "Solves a time-dependent linear PDE with discontinuous right-hand side.\n";
2fa9584fbSIlya Fursov
3fa9584fbSIlya Fursov /* ------------------------------------------------------------------------
4fa9584fbSIlya Fursov
5fa9584fbSIlya Fursov This program solves the one-dimensional quench front problem modeling a cooled
6fa9584fbSIlya Fursov liquid rising on a hot metal rod
7fa9584fbSIlya Fursov u_t = u_xx + g(u),
8fa9584fbSIlya Fursov with
9fa9584fbSIlya Fursov g(u) = -Au if u <= u_c,
10fa9584fbSIlya Fursov = 0 if u > u_c
11fa9584fbSIlya Fursov on the domain 0 <= x <= 1, with the boundary conditions
12fa9584fbSIlya Fursov u(t,0) = 0, u_x(t,1) = 0,
13fa9584fbSIlya Fursov and the initial condition
14fa9584fbSIlya Fursov u(0,x) = 0 if 0 <= x <= 0.1,
15fa9584fbSIlya Fursov = (x - 0.1)/0.15 if 0.1 < x < 0.25
16fa9584fbSIlya Fursov = 1 if 0.25 <= x <= 1
17fa9584fbSIlya Fursov We discretize the right-hand side using finite differences with
18fa9584fbSIlya Fursov uniform grid spacing h:
19fa9584fbSIlya Fursov u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
20fa9584fbSIlya Fursov
21fa9584fbSIlya Fursov Reference: L. Shampine and S. Thompson, "Event Location for Ordinary Differential Equations",
22fa9584fbSIlya Fursov http://www.radford.edu/~thompson/webddes/eventsweb.pdf
23fa9584fbSIlya Fursov ------------------------------------------------------------------------- */
24fa9584fbSIlya Fursov
25fa9584fbSIlya Fursov #include <petscdmda.h>
26fa9584fbSIlya Fursov #include <petscts.h>
27fa9584fbSIlya Fursov /*
28fa9584fbSIlya Fursov User-defined application context - contains data needed by the
29fa9584fbSIlya Fursov application-provided call-back routines.
30fa9584fbSIlya Fursov */
31fa9584fbSIlya Fursov typedef struct {
32fa9584fbSIlya Fursov PetscReal A;
33fa9584fbSIlya Fursov PetscReal uc;
34fa9584fbSIlya Fursov PetscInt *sw;
35fa9584fbSIlya Fursov } AppCtx;
36fa9584fbSIlya Fursov
InitialConditions(Vec U,DM da,AppCtx * app)37fa9584fbSIlya Fursov PetscErrorCode InitialConditions(Vec U, DM da, AppCtx *app)
38fa9584fbSIlya Fursov {
39fa9584fbSIlya Fursov Vec xcoord;
40fa9584fbSIlya Fursov PetscScalar *x, *u;
41fa9584fbSIlya Fursov PetscInt lsize, M, xs, xm, i;
42fa9584fbSIlya Fursov
43fa9584fbSIlya Fursov PetscFunctionBeginUser;
44fa9584fbSIlya Fursov PetscCall(DMGetCoordinates(da, &xcoord));
45fa9584fbSIlya Fursov PetscCall(DMDAVecGetArrayRead(da, xcoord, &x));
46fa9584fbSIlya Fursov
47fa9584fbSIlya Fursov PetscCall(VecGetLocalSize(U, &lsize));
48fa9584fbSIlya Fursov PetscCall(PetscMalloc1(lsize, &app->sw));
49fa9584fbSIlya Fursov
50fa9584fbSIlya Fursov PetscCall(DMDAVecGetArray(da, U, &u));
51fa9584fbSIlya Fursov
52fa9584fbSIlya Fursov PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
53fa9584fbSIlya Fursov PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
54fa9584fbSIlya Fursov
55fa9584fbSIlya Fursov for (i = xs; i < xs + xm; i++) {
56fa9584fbSIlya Fursov if (x[i] <= 0.1) u[i] = 0.;
57fa9584fbSIlya Fursov else if (x[i] > 0.1 && x[i] < 0.25) u[i] = (x[i] - 0.1) / 0.15;
58fa9584fbSIlya Fursov else u[i] = 1.0;
59fa9584fbSIlya Fursov
60fa9584fbSIlya Fursov app->sw[i - xs] = 1;
61fa9584fbSIlya Fursov }
62fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArray(da, U, &u));
63fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArrayRead(da, xcoord, &x));
64fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
65fa9584fbSIlya Fursov }
66fa9584fbSIlya Fursov
EventFunction(TS ts,PetscReal t,Vec U,PetscReal * fvalue,PetscCtx ctx)67*2a8381b2SBarry Smith PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, PetscCtx ctx)
68fa9584fbSIlya Fursov {
69fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx;
70fa9584fbSIlya Fursov const PetscScalar *u;
71fa9584fbSIlya Fursov PetscInt i, lsize;
72fa9584fbSIlya Fursov
73fa9584fbSIlya Fursov PetscFunctionBeginUser;
74fa9584fbSIlya Fursov PetscCall(VecGetLocalSize(U, &lsize));
75fa9584fbSIlya Fursov PetscCall(VecGetArrayRead(U, &u));
76fa9584fbSIlya Fursov for (i = 0; i < lsize; i++) fvalue[i] = PetscRealPart(u[i]) - app->uc;
77fa9584fbSIlya Fursov PetscCall(VecRestoreArrayRead(U, &u));
78fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
79fa9584fbSIlya Fursov }
80fa9584fbSIlya Fursov
PostEventFunction(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,PetscCtx ctx)81*2a8381b2SBarry Smith PetscErrorCode PostEventFunction(TS ts, PetscInt nevents_zero, PetscInt events_zero[], PetscReal t, Vec U, PetscBool forwardsolve, PetscCtx ctx)
82fa9584fbSIlya Fursov {
83fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx;
84fa9584fbSIlya Fursov PetscInt i, idx;
85fa9584fbSIlya Fursov
86fa9584fbSIlya Fursov PetscFunctionBeginUser;
87fa9584fbSIlya Fursov for (i = 0; i < nevents_zero; i++) {
88fa9584fbSIlya Fursov idx = events_zero[i];
89fa9584fbSIlya Fursov app->sw[idx] = 0;
90fa9584fbSIlya Fursov }
91fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
92fa9584fbSIlya Fursov }
93fa9584fbSIlya Fursov
94fa9584fbSIlya Fursov /*
95fa9584fbSIlya Fursov Defines the ODE passed to the ODE solver
96fa9584fbSIlya Fursov */
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)97*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
98fa9584fbSIlya Fursov {
99fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx;
100fa9584fbSIlya Fursov PetscScalar *f;
101fa9584fbSIlya Fursov const PetscScalar *u, *udot;
102fa9584fbSIlya Fursov DM da;
103fa9584fbSIlya Fursov PetscInt M, xs, xm, i;
104fa9584fbSIlya Fursov PetscReal h, h2;
105fa9584fbSIlya Fursov Vec Ulocal;
106fa9584fbSIlya Fursov
107fa9584fbSIlya Fursov PetscFunctionBeginUser;
108fa9584fbSIlya Fursov PetscCall(TSGetDM(ts, &da));
109fa9584fbSIlya Fursov
110fa9584fbSIlya Fursov PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
111fa9584fbSIlya Fursov PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
112fa9584fbSIlya Fursov
113fa9584fbSIlya Fursov PetscCall(DMGetLocalVector(da, &Ulocal));
114fa9584fbSIlya Fursov PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, Ulocal));
115fa9584fbSIlya Fursov PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, Ulocal));
116fa9584fbSIlya Fursov
117fa9584fbSIlya Fursov h = 1.0 / (M - 1);
118fa9584fbSIlya Fursov h2 = h * h;
119fa9584fbSIlya Fursov PetscCall(DMDAVecGetArrayRead(da, Udot, &udot));
120fa9584fbSIlya Fursov PetscCall(DMDAVecGetArrayRead(da, Ulocal, &u));
121fa9584fbSIlya Fursov PetscCall(DMDAVecGetArray(da, F, &f));
122fa9584fbSIlya Fursov
123fa9584fbSIlya Fursov for (i = xs; i < xs + xm; i++) {
124fa9584fbSIlya Fursov if (i == 0) {
125fa9584fbSIlya Fursov f[i] = u[i];
126fa9584fbSIlya Fursov } else if (i == M - 1) {
127fa9584fbSIlya Fursov f[i] = (u[i] - u[i - 1]) / h;
128fa9584fbSIlya Fursov } else {
129fa9584fbSIlya Fursov f[i] = (u[i + 1] - 2 * u[i] + u[i - 1]) / h2 + app->sw[i - xs] * (-app->A * u[i]) - udot[i];
130fa9584fbSIlya Fursov }
131fa9584fbSIlya Fursov }
132fa9584fbSIlya Fursov
133fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot));
134fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArrayRead(da, Ulocal, &u));
135fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArray(da, F, &f));
136fa9584fbSIlya Fursov PetscCall(DMRestoreLocalVector(da, &Ulocal));
137fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
138fa9584fbSIlya Fursov }
139fa9584fbSIlya Fursov
140fa9584fbSIlya Fursov /*
141fa9584fbSIlya Fursov Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
142fa9584fbSIlya Fursov */
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)143*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
144fa9584fbSIlya Fursov {
145fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx;
146fa9584fbSIlya Fursov DM da;
147fa9584fbSIlya Fursov MatStencil row, col[3];
148fa9584fbSIlya Fursov PetscScalar v[3];
149fa9584fbSIlya Fursov PetscInt M, xs, xm, i;
150fa9584fbSIlya Fursov PetscReal h, h2;
151fa9584fbSIlya Fursov
152fa9584fbSIlya Fursov PetscFunctionBeginUser;
153fa9584fbSIlya Fursov PetscCall(TSGetDM(ts, &da));
154fa9584fbSIlya Fursov
155fa9584fbSIlya Fursov PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
156fa9584fbSIlya Fursov PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
157fa9584fbSIlya Fursov
158fa9584fbSIlya Fursov h = 1.0 / (M - 1);
159fa9584fbSIlya Fursov h2 = h * h;
160fa9584fbSIlya Fursov for (i = xs; i < xs + xm; i++) {
161fa9584fbSIlya Fursov row.i = i;
162fa9584fbSIlya Fursov if (i == 0) {
163fa9584fbSIlya Fursov v[0] = 1.0;
164fa9584fbSIlya Fursov PetscCall(MatSetValuesStencil(A, 1, &row, 1, &row, v, INSERT_VALUES));
165fa9584fbSIlya Fursov } else if (i == M - 1) {
166fa9584fbSIlya Fursov col[0].i = i;
167fa9584fbSIlya Fursov v[0] = 1 / h;
168fa9584fbSIlya Fursov col[1].i = i - 1;
169fa9584fbSIlya Fursov v[1] = -1 / h;
170fa9584fbSIlya Fursov PetscCall(MatSetValuesStencil(A, 1, &row, 2, col, v, INSERT_VALUES));
171fa9584fbSIlya Fursov } else {
172fa9584fbSIlya Fursov col[0].i = i + 1;
173fa9584fbSIlya Fursov v[0] = 1 / h2;
174fa9584fbSIlya Fursov col[1].i = i;
175fa9584fbSIlya Fursov v[1] = -2 / h2 + app->sw[i - xs] * (-app->A) - a;
176fa9584fbSIlya Fursov col[2].i = i - 1;
177fa9584fbSIlya Fursov v[2] = 1 / h2;
178fa9584fbSIlya Fursov PetscCall(MatSetValuesStencil(A, 1, &row, 3, col, v, INSERT_VALUES));
179fa9584fbSIlya Fursov }
180fa9584fbSIlya Fursov }
181fa9584fbSIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
182fa9584fbSIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
183fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
184fa9584fbSIlya Fursov }
185fa9584fbSIlya Fursov
main(int argc,char ** argv)186fa9584fbSIlya Fursov int main(int argc, char **argv)
187fa9584fbSIlya Fursov {
188fa9584fbSIlya Fursov TS ts; /* ODE integrator */
189fa9584fbSIlya Fursov Vec U; /* solution will be stored here */
190fa9584fbSIlya Fursov Mat J; /* Jacobian matrix */
191fa9584fbSIlya Fursov PetscInt n = 16;
192fa9584fbSIlya Fursov AppCtx app;
193fa9584fbSIlya Fursov DM da;
194fa9584fbSIlya Fursov
195fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196fa9584fbSIlya Fursov Initialize program
197fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198fa9584fbSIlya Fursov PetscFunctionBeginUser;
199c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
200fa9584fbSIlya Fursov
201fa9584fbSIlya Fursov PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex22 options", "");
202fa9584fbSIlya Fursov {
203fa9584fbSIlya Fursov app.A = 200000;
204fa9584fbSIlya Fursov PetscCall(PetscOptionsReal("-A", "", "", app.A, &app.A, NULL));
205fa9584fbSIlya Fursov app.uc = 0.5;
206fa9584fbSIlya Fursov PetscCall(PetscOptionsReal("-uc", "", "", app.uc, &app.uc, NULL));
207fa9584fbSIlya Fursov }
208fa9584fbSIlya Fursov PetscOptionsEnd();
209fa9584fbSIlya Fursov
210fa9584fbSIlya Fursov PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, -n, 1, 1, 0, &da));
211fa9584fbSIlya Fursov PetscCall(DMSetFromOptions(da));
212fa9584fbSIlya Fursov PetscCall(DMSetUp(da));
213fa9584fbSIlya Fursov PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0, 0, 0, 0));
214fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215fa9584fbSIlya Fursov Create necessary matrix and vectors
216fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217fa9584fbSIlya Fursov PetscCall(DMCreateMatrix(da, &J));
218fa9584fbSIlya Fursov PetscCall(DMCreateGlobalVector(da, &U));
219fa9584fbSIlya Fursov
220fa9584fbSIlya Fursov PetscCall(InitialConditions(U, da, &app));
221fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222fa9584fbSIlya Fursov Create timestepping solver context
223fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224fa9584fbSIlya Fursov PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
225fa9584fbSIlya Fursov PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
226fa9584fbSIlya Fursov PetscCall(TSSetType(ts, TSROSW));
2278434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, (void *)&app));
2288434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobianFn *)IJacobian, (void *)&app));
229fa9584fbSIlya Fursov
230fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231fa9584fbSIlya Fursov Set initial conditions
232fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233fa9584fbSIlya Fursov PetscCall(TSSetSolution(ts, U));
234fa9584fbSIlya Fursov
235fa9584fbSIlya Fursov PetscCall(TSSetDM(ts, da));
236fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237fa9584fbSIlya Fursov Set solver options
238fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239fa9584fbSIlya Fursov PetscCall(TSSetTimeStep(ts, 0.1));
240fa9584fbSIlya Fursov PetscCall(TSSetMaxTime(ts, 30.0));
241fa9584fbSIlya Fursov PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
242fa9584fbSIlya Fursov PetscCall(TSSetFromOptions(ts));
243fa9584fbSIlya Fursov
244fa9584fbSIlya Fursov PetscInt lsize;
245fa9584fbSIlya Fursov PetscCall(VecGetLocalSize(U, &lsize));
246fa9584fbSIlya Fursov PetscInt *direction;
247fa9584fbSIlya Fursov PetscBool *terminate;
248fa9584fbSIlya Fursov PetscInt i;
249fa9584fbSIlya Fursov PetscCall(PetscMalloc1(lsize, &direction));
250fa9584fbSIlya Fursov PetscCall(PetscMalloc1(lsize, &terminate));
251fa9584fbSIlya Fursov for (i = 0; i < lsize; i++) {
252fa9584fbSIlya Fursov direction[i] = -1;
253fa9584fbSIlya Fursov terminate[i] = PETSC_FALSE;
254fa9584fbSIlya Fursov }
255fa9584fbSIlya Fursov PetscCall(TSSetEventHandler(ts, lsize, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
256fa9584fbSIlya Fursov
257fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258fa9584fbSIlya Fursov Run timestepping solver
259fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260fa9584fbSIlya Fursov PetscCall(TSSolve(ts, U));
261fa9584fbSIlya Fursov
262fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263fa9584fbSIlya Fursov Free work space. All PETSc objects should be destroyed when they are no longer needed.
264fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
265fa9584fbSIlya Fursov
266fa9584fbSIlya Fursov PetscCall(MatDestroy(&J));
267fa9584fbSIlya Fursov PetscCall(VecDestroy(&U));
268fa9584fbSIlya Fursov PetscCall(DMDestroy(&da));
269fa9584fbSIlya Fursov PetscCall(TSDestroy(&ts));
270fa9584fbSIlya Fursov PetscCall(PetscFree(direction));
271fa9584fbSIlya Fursov PetscCall(PetscFree(terminate));
272fa9584fbSIlya Fursov
273fa9584fbSIlya Fursov PetscCall(PetscFree(app.sw));
274fa9584fbSIlya Fursov PetscCall(PetscFinalize());
275fa9584fbSIlya Fursov return 0;
276fa9584fbSIlya Fursov }
277