xref: /honee/src/setupts.c (revision bb8a0c61f21224cefcdd60e71004bb99df1e9a58)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Time-stepping functions for Navier-Stokes example using PETSc
10 
11 #include "../navierstokes.h"
12 #include "../qfunctions/mass.h"
13 
14 // Compute mass matrix for explicit scheme
15 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data,
16                                        Vec M) {
17   Vec            M_loc;
18   CeedQFunction  qf_mass;
19   CeedOperator   op_mass;
20   CeedVector     m_ceed, ones_vec;
21   CeedInt        num_comp_q, q_data_size;
22   PetscErrorCode ierr;
23   PetscFunctionBeginUser;
24 
25   // CEED Restriction
26   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
27   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
28   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
29   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
30   CeedVectorSetValue(ones_vec, 1.0);
31 
32   // CEED QFunction
33   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
34   CeedQFunctionAddInput(qf_mass, "q", num_comp_q, CEED_EVAL_INTERP);
35   CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE);
36   CeedQFunctionAddOutput(qf_mass, "v", num_comp_q, CEED_EVAL_INTERP);
37 
38   // CEED Operator
39   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
40   CeedOperatorSetField(op_mass, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
41                        CEED_VECTOR_ACTIVE);
42   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i,
43                        CEED_BASIS_COLLOCATED, ceed_data->q_data);
44   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
45                        CEED_VECTOR_ACTIVE);
46 
47   // Place PETSc vector in CEED vector
48   CeedScalar *m;
49   PetscMemType m_mem_type;
50   ierr = DMGetLocalVector(dm, &M_loc); CHKERRQ(ierr);
51   ierr = VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type);
52   CHKERRQ(ierr);
53   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
54 
55   // Apply CEED Operator
56   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
57 
58   // Restore vectors
59   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
60   ierr = VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m);
61   CHKERRQ(ierr);
62 
63   // Local-to-Global
64   ierr = VecZeroEntries(M); CHKERRQ(ierr);
65   ierr = DMLocalToGlobal(dm, M_loc, ADD_VALUES, M); CHKERRQ(ierr);
66   ierr = DMRestoreLocalVector(dm, &M_loc); CHKERRQ(ierr);
67 
68   // Invert diagonally lumped mass vector for RHS function
69   ierr = VecReciprocal(M); CHKERRQ(ierr);
70 
71   // Cleanup
72   CeedVectorDestroy(&ones_vec);
73   CeedVectorDestroy(&m_ceed);
74   CeedQFunctionDestroy(&qf_mass);
75   CeedOperatorDestroy(&op_mass);
76 
77   PetscFunctionReturn(0);
78 }
79 
80 // RHS (Explicit time-stepper) function setup
81 //   This is the RHS of the ODE, given as u_t = G(t,u)
82 //   This function takes in a state vector Q and writes into G
83 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
84 
85   User           user = *(User *)user_data;
86   PetscScalar    *q, *g;
87   Vec            Q_loc, G_loc;
88   PetscMemType   q_mem_type, g_mem_type;
89   PetscErrorCode ierr;
90   PetscFunctionBeginUser;
91 
92   // Update context field labels
93   if (user->phys->solution_time_label)
94     CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t);
95   if (user->phys->timestep_size_label) {
96     PetscScalar dt;
97     ierr = TSGetTimeStep(ts,&dt); CHKERRQ(ierr);
98     CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label,
99                                  &dt);
100   }
101 
102   // Get local vectors
103   ierr = DMGetLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
104   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
105 
106   // Global-to-local
107   ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr);
108   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
109   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, 0.0,
110                                     NULL, NULL, NULL); CHKERRQ(ierr);
111   ierr = VecZeroEntries(G_loc); CHKERRQ(ierr);
112 
113   // Place PETSc vectors in CEED vectors
114   ierr = VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type);
115   CHKERRQ(ierr);
116   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
117   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
118   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
119 
120   // Apply CEED operator
121   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed,
122                     CEED_REQUEST_IMMEDIATE);
123 
124   // Restore vectors
125   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
126   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
127   ierr = VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q);
128   CHKERRQ(ierr);
129   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
130 
131   // Local-to-Global
132   ierr = VecZeroEntries(G); CHKERRQ(ierr);
133   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
134 
135   // Inverse of the lumped mass matrix (M is Minv)
136   ierr = VecPointwiseMult(G, G, user->M); CHKERRQ(ierr);
137 
138   // Restore vectors
139   ierr = DMRestoreLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
140   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
141 
142   PetscFunctionReturn(0);
143 }
144 
145 // Implicit time-stepper function setup
146 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G,
147                             void *user_data) {
148   User              user = *(User *)user_data;
149   const PetscScalar *q, *q_dot;
150   PetscScalar       *g;
151   Vec               Q_loc, Q_dot_loc, G_loc;
152   PetscMemType      q_mem_type, q_dot_mem_type, g_mem_type;
153   PetscErrorCode    ierr;
154   PetscFunctionBeginUser;
155 
156   // Update context field labels
157   if (user->phys->solution_time_label)
158     CeedOperatorContextSetDouble(user->op_ifunction,
159                                  user->phys->solution_time_label, &t);
160   if (user->phys->timestep_size_label) {
161     PetscScalar dt;
162     ierr = TSGetTimeStep(ts,&dt); CHKERRQ(ierr);
163     CeedOperatorContextSetDouble(user->op_ifunction,
164                                  user->phys->timestep_size_label, &dt);
165   }
166 
167   // Get local vectors
168   ierr = DMGetLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
169   ierr = DMGetLocalVector(user->dm, &Q_dot_loc); CHKERRQ(ierr);
170   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
171 
172   // Global-to-local
173   ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr);
174   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
175   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, 0.0,
176                                     NULL, NULL, NULL); CHKERRQ(ierr);
177   ierr = VecZeroEntries(Q_dot_loc); CHKERRQ(ierr);
178   ierr = DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc);
179   CHKERRQ(ierr);
180   ierr = VecZeroEntries(G_loc); CHKERRQ(ierr);
181 
182   // Place PETSc vectors in CEED vectors
183   ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr);
184   ierr = VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type);
185   CHKERRQ(ierr);
186   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
187   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER,
188                      (PetscScalar *)q);
189   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type),
190                      CEED_USE_POINTER, (PetscScalar *)q_dot);
191   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
192 
193   // Apply CEED operator
194   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed,
195                     CEED_REQUEST_IMMEDIATE);
196 
197   // Restore vectors
198   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
199   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
200   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
201   ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr);
202   ierr = VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot); CHKERRQ(ierr);
203   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
204 
205   // Local-to-Global
206   ierr = VecZeroEntries(G); CHKERRQ(ierr);
207   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
208 
209   // Restore vectors
210   ierr = DMRestoreLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
211   ierr = DMRestoreLocalVector(user->dm, &Q_dot_loc); CHKERRQ(ierr);
212   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
213 
214   PetscFunctionReturn(0);
215 }
216 
217 // User provided TS Monitor
218 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time,
219                             Vec Q, void *ctx) {
220   User           user = ctx;
221   Vec            Q_loc;
222   char           file_path[PETSC_MAX_PATH_LEN];
223   PetscViewer    viewer;
224   PetscErrorCode ierr;
225   PetscFunctionBeginUser;
226 
227   // Print every 'output_freq' steps
228   if (step_no % user->app_ctx->output_freq != 0)
229     PetscFunctionReturn(0);
230 
231   // Set up output
232   ierr = DMGetLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
233   ierr = PetscObjectSetName((PetscObject)Q_loc, "StateVec"); CHKERRQ(ierr);
234   ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr);
235   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
236 
237   // Output
238   ierr = PetscSNPrintf(file_path, sizeof file_path,
239                        "%s/ns-%03" PetscInt_FMT ".vtu",
240                        user->app_ctx->output_dir, step_no + user->app_ctx->cont_steps);
241   CHKERRQ(ierr);
242   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path,
243                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
244   ierr = VecView(Q_loc, viewer); CHKERRQ(ierr);
245   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
246   if (user->dm_viz) {
247     Vec         Q_refined, Q_refined_loc;
248     char        file_path_refined[PETSC_MAX_PATH_LEN];
249     PetscViewer viewer_refined;
250 
251     ierr = DMGetGlobalVector(user->dm_viz, &Q_refined); CHKERRQ(ierr);
252     ierr = DMGetLocalVector(user->dm_viz, &Q_refined_loc); CHKERRQ(ierr);
253     ierr = PetscObjectSetName((PetscObject)Q_refined_loc, "Refined");
254     CHKERRQ(ierr);
255     ierr = MatInterpolate(user->interp_viz, Q, Q_refined); CHKERRQ(ierr);
256     ierr = VecZeroEntries(Q_refined_loc); CHKERRQ(ierr);
257     ierr = DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc);
258     CHKERRQ(ierr);
259     ierr = PetscSNPrintf(file_path_refined, sizeof file_path_refined,
260                          "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir,
261                          step_no + user->app_ctx->cont_steps);
262     CHKERRQ(ierr);
263     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined),
264                               file_path_refined, FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
265     ierr = VecView(Q_refined_loc, viewer_refined); CHKERRQ(ierr);
266     ierr = DMRestoreLocalVector(user->dm_viz, &Q_refined_loc); CHKERRQ(ierr);
267     ierr = DMRestoreGlobalVector(user->dm_viz, &Q_refined); CHKERRQ(ierr);
268     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
269   }
270   ierr = DMRestoreLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
271 
272   // Save data in a binary file for continuation of simulations
273   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin",
274                        user->app_ctx->output_dir); CHKERRQ(ierr);
275   ierr = PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer);
276   CHKERRQ(ierr);
277   ierr = VecView(Q, viewer); CHKERRQ(ierr);
278   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
279 
280   // Save time stamp
281   // Dimensionalize time back
282   time /= user->units->second;
283   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin",
284                        user->app_ctx->output_dir); CHKERRQ(ierr);
285   ierr = PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer);
286   CHKERRQ(ierr);
287   #if PETSC_VERSION_GE(3,13,0)
288   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
289   #else
290   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
291   #endif
292   CHKERRQ(ierr);
293   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
294 
295   PetscFunctionReturn(0);
296 }
297 
298 // TS: Create, setup, and solve
299 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys,
300                           Vec *Q, PetscScalar *f_time, TS *ts) {
301   MPI_Comm       comm = user->comm;
302   TSAdapt        adapt;
303   PetscScalar    final_time;
304   PetscErrorCode ierr;
305   PetscFunctionBeginUser;
306 
307   ierr = TSCreate(comm, ts); CHKERRQ(ierr);
308   ierr = TSSetDM(*ts, dm); CHKERRQ(ierr);
309   if (phys->implicit) {
310     ierr = TSSetType(*ts, TSBDF); CHKERRQ(ierr);
311     if (user->op_ifunction) {
312       ierr = TSSetIFunction(*ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
313     } else { // Implicit integrators can fall back to using an RHSFunction
314       ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
315     }
316   } else {
317     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
318                                  "Problem does not provide RHSFunction");
319     ierr = TSSetType(*ts, TSRK); CHKERRQ(ierr);
320     ierr = TSRKSetType(*ts, TSRK5F); CHKERRQ(ierr);
321     ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
322   }
323   ierr = TSSetMaxTime(*ts, 500. * user->units->second); CHKERRQ(ierr);
324   ierr = TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
325   ierr = TSSetTimeStep(*ts, 1.e-2 * user->units->second); CHKERRQ(ierr);
326   if (app_ctx->test_mode) {ierr = TSSetMaxSteps(*ts, 10); CHKERRQ(ierr);}
327   ierr = TSGetAdapt(*ts, &adapt); CHKERRQ(ierr);
328   ierr = TSAdaptSetStepLimits(adapt,
329                               1.e-12 * user->units->second,
330                               1.e2 * user->units->second); CHKERRQ(ierr);
331   ierr = TSSetFromOptions(*ts); CHKERRQ(ierr);
332   if (!app_ctx->cont_steps) { // print initial condition
333     if (!app_ctx->test_mode) {
334       ierr = TSMonitor_NS(*ts, 0, 0., *Q, user); CHKERRQ(ierr);
335     }
336   } else { // continue from time of last output
337     PetscReal   time;
338     PetscInt    count;
339     PetscViewer viewer;
340     char        file_path[PETSC_MAX_PATH_LEN];
341     ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin",
342                          app_ctx->output_dir); CHKERRQ(ierr);
343     ierr = PetscViewerBinaryOpen(comm, file_path, FILE_MODE_READ, &viewer);
344     CHKERRQ(ierr);
345     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
346     CHKERRQ(ierr);
347     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
348     ierr = TSSetTime(*ts, time * user->units->second); CHKERRQ(ierr);
349   }
350   if (!app_ctx->test_mode) {
351     ierr = TSMonitorSet(*ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
352   }
353 
354   // Solve
355   double start, cpu_time_used;
356   start = MPI_Wtime();
357   ierr = PetscBarrier((PetscObject) *ts); CHKERRQ(ierr);
358   ierr = TSSolve(*ts, *Q); CHKERRQ(ierr);
359   cpu_time_used = MPI_Wtime() - start;
360   ierr = TSGetSolveTime(*ts, &final_time); CHKERRQ(ierr);
361   *f_time = final_time;
362   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
363                        comm); CHKERRQ(ierr);
364   if (!app_ctx->test_mode) {
365     ierr = PetscPrintf(PETSC_COMM_WORLD,
366                        "Time taken for solution (sec): %g\n",
367                        (double)cpu_time_used); CHKERRQ(ierr);
368   }
369   PetscFunctionReturn(0);
370 }
371