xref: /honee/src/setupts.c (revision 9f59f36e7a2e19e53ad001515b81028b86a73a2c)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11a515125bSLeila Ghaffari #include "../navierstokes.h"
12a515125bSLeila Ghaffari 
13a515125bSLeila Ghaffari // Compute mass matrix for explicit scheme
142b916ea7SJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
15a515125bSLeila Ghaffari   Vec           M_loc;
16a515125bSLeila Ghaffari   CeedQFunction qf_mass;
17a515125bSLeila Ghaffari   CeedOperator  op_mass;
18a515125bSLeila Ghaffari   CeedVector    m_ceed, ones_vec;
19a515125bSLeila Ghaffari   CeedInt       num_comp_q, q_data_size;
20a515125bSLeila Ghaffari   PetscFunctionBeginUser;
21a515125bSLeila Ghaffari 
22a515125bSLeila Ghaffari   // CEED Restriction
23a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
24a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
25a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
26a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
27a515125bSLeila Ghaffari   CeedVectorSetValue(ones_vec, 1.0);
28a515125bSLeila Ghaffari 
29a515125bSLeila Ghaffari   // CEED QFunction
30*9f59f36eSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
31a515125bSLeila Ghaffari 
32a515125bSLeila Ghaffari   // CEED Operator
33a515125bSLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
34*9f59f36eSJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
352b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
362b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
37a515125bSLeila Ghaffari 
38a515125bSLeila Ghaffari   // Place PETSc vector in CEED vector
39a515125bSLeila Ghaffari   CeedScalar  *m;
40a515125bSLeila Ghaffari   PetscMemType m_mem_type;
412b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &M_loc));
422b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type));
43a515125bSLeila Ghaffari   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
44a515125bSLeila Ghaffari 
45a515125bSLeila Ghaffari   // Apply CEED Operator
46a515125bSLeila Ghaffari   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
47a515125bSLeila Ghaffari 
48a515125bSLeila Ghaffari   // Restore vectors
49a515125bSLeila Ghaffari   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
502b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m));
51a515125bSLeila Ghaffari 
52a515125bSLeila Ghaffari   // Local-to-Global
532b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(M));
542b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M));
552b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &M_loc));
56a515125bSLeila Ghaffari 
57a515125bSLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
582b916ea7SJeremy L Thompson   PetscCall(VecReciprocal(M));
59a515125bSLeila Ghaffari 
60a515125bSLeila Ghaffari   // Cleanup
61a515125bSLeila Ghaffari   CeedVectorDestroy(&ones_vec);
62a515125bSLeila Ghaffari   CeedVectorDestroy(&m_ceed);
63a515125bSLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
64a515125bSLeila Ghaffari   CeedOperatorDestroy(&op_mass);
65a515125bSLeila Ghaffari 
66a515125bSLeila Ghaffari   PetscFunctionReturn(0);
67a515125bSLeila Ghaffari }
68a515125bSLeila Ghaffari 
69a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
70a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
71a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
72a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
73a515125bSLeila Ghaffari   User         user = *(User *)user_data;
74a515125bSLeila Ghaffari   PetscScalar *q, *g;
75e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
76a515125bSLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
77a515125bSLeila Ghaffari   PetscFunctionBeginUser;
78a515125bSLeila Ghaffari 
79e2f84137SJeremy L Thompson   // Get local vector
802b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
81e2f84137SJeremy L Thompson 
82e2f84137SJeremy L Thompson   // Update time dependent data
83e2f84137SJeremy L Thompson   if (user->time != t) {
842b916ea7SJeremy L Thompson     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
85e2f84137SJeremy L Thompson     if (user->phys->solution_time_label) {
8692ada588SJames Wright       CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t);
87e2f84137SJeremy L Thompson     }
88e2f84137SJeremy L Thompson     user->time = t;
89e2f84137SJeremy L Thompson   }
90bb8a0c61SJames Wright   if (user->phys->timestep_size_label) {
91bb8a0c61SJames Wright     PetscScalar dt;
922b916ea7SJeremy L Thompson     PetscCall(TSGetTimeStep(ts, &dt));
93e2f84137SJeremy L Thompson     if (user->dt != dt) {
942b916ea7SJeremy L Thompson       CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label, &dt);
95e2f84137SJeremy L Thompson       user->dt = dt;
96e2f84137SJeremy L Thompson     }
97bb8a0c61SJames Wright   }
98a515125bSLeila Ghaffari 
99a515125bSLeila Ghaffari   // Global-to-local
1002b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
101a515125bSLeila Ghaffari 
102a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
1032b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type));
1042b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
105a515125bSLeila Ghaffari   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
106a515125bSLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
107a515125bSLeila Ghaffari 
108a515125bSLeila Ghaffari   // Apply CEED operator
1092b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
110a515125bSLeila Ghaffari 
111a515125bSLeila Ghaffari   // Restore vectors
112a515125bSLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
113a515125bSLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
1142b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q));
1152b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
116a515125bSLeila Ghaffari 
117a515125bSLeila Ghaffari   // Local-to-Global
1182b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1192b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
120a515125bSLeila Ghaffari 
121a515125bSLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1222b916ea7SJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
123a515125bSLeila Ghaffari 
124a515125bSLeila Ghaffari   // Restore vectors
1252b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
126a515125bSLeila Ghaffari 
127a515125bSLeila Ghaffari   PetscFunctionReturn(0);
128a515125bSLeila Ghaffari }
129a515125bSLeila Ghaffari 
130a515125bSLeila Ghaffari // Implicit time-stepper function setup
1312b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
132a515125bSLeila Ghaffari   User               user = *(User *)user_data;
133a515125bSLeila Ghaffari   const PetscScalar *q, *q_dot;
134a515125bSLeila Ghaffari   PetscScalar       *g;
135e2f84137SJeremy L Thompson   Vec                Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
136a515125bSLeila Ghaffari   PetscMemType       q_mem_type, q_dot_mem_type, g_mem_type;
137a515125bSLeila Ghaffari   PetscFunctionBeginUser;
138a515125bSLeila Ghaffari 
139e2f84137SJeremy L Thompson   // Get local vectors
1402b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
141e2f84137SJeremy L Thompson 
142e2f84137SJeremy L Thompson   // Update time dependent data
143e2f84137SJeremy L Thompson   if (user->time != t) {
1442b916ea7SJeremy L Thompson     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
145e2f84137SJeremy L Thompson     if (user->phys->solution_time_label) {
1462b916ea7SJeremy L Thompson       CeedOperatorContextSetDouble(user->op_ifunction, user->phys->solution_time_label, &t);
147e2f84137SJeremy L Thompson     }
148e2f84137SJeremy L Thompson     user->time = t;
149e2f84137SJeremy L Thompson   }
150bb8a0c61SJames Wright   if (user->phys->timestep_size_label) {
151bb8a0c61SJames Wright     PetscScalar dt;
1522b916ea7SJeremy L Thompson     PetscCall(TSGetTimeStep(ts, &dt));
153e2f84137SJeremy L Thompson     if (user->dt != dt) {
1542b916ea7SJeremy L Thompson       CeedOperatorContextSetDouble(user->op_ifunction, user->phys->timestep_size_label, &dt);
155e2f84137SJeremy L Thompson       user->dt = dt;
156e2f84137SJeremy L Thompson     }
157bb8a0c61SJames Wright   }
158a515125bSLeila Ghaffari 
159a515125bSLeila Ghaffari   // Global-to-local
1602b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
1612b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
162a515125bSLeila Ghaffari 
163a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
1642b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
1652b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type));
1662b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
1672b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
1682b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), CEED_USE_POINTER, (PetscScalar *)q_dot);
169a515125bSLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
170a515125bSLeila Ghaffari 
171a515125bSLeila Ghaffari   // Apply CEED operator
1722b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
173a515125bSLeila Ghaffari 
174a515125bSLeila Ghaffari   // Restore vectors
175a515125bSLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
176a515125bSLeila Ghaffari   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
177a515125bSLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
1782b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
1792b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot));
1802b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
181a515125bSLeila Ghaffari 
182a515125bSLeila Ghaffari   // Local-to-Global
1832b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1842b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
185a515125bSLeila Ghaffari 
186a515125bSLeila Ghaffari   // Restore vectors
1872b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
188a515125bSLeila Ghaffari 
189a515125bSLeila Ghaffari   PetscFunctionReturn(0);
190a515125bSLeila Ghaffari }
191a515125bSLeila Ghaffari 
192f0b65372SJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
193f0b65372SJed Brown   User               user;
194f0b65372SJed Brown   const PetscScalar *q;
195f0b65372SJed Brown   PetscScalar       *g;
196f0b65372SJed Brown   PetscMemType       q_mem_type, g_mem_type;
197f0b65372SJed Brown   PetscFunctionBeginUser;
1982b916ea7SJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
199e2f84137SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
200e2f84137SJeremy L Thompson       G_loc;
201e2f84137SJeremy L Thompson 
202f0b65372SJed Brown   // Get local vectors
2032b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
204f0b65372SJed Brown 
205f0b65372SJed Brown   // Global-to-local
2062b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
207f0b65372SJed Brown 
208f0b65372SJed Brown   // Place PETSc vectors in CEED vectors
2092b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2102b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2112b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
212f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
213f0b65372SJed Brown 
214f0b65372SJed Brown   // Apply CEED operator
2152b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
216f0b65372SJed Brown 
217f0b65372SJed Brown   // Restore vectors
218f0b65372SJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
219f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2202b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2212b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
222f0b65372SJed Brown 
223f0b65372SJed Brown   // Local-to-Global
2242b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2252b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
226f0b65372SJed Brown 
227f0b65372SJed Brown   // Restore vectors
2282b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
229f0b65372SJed Brown   PetscFunctionReturn(0);
230f0b65372SJed Brown }
231f0b65372SJed Brown 
232f0b65372SJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
233f0b65372SJed Brown   User         user;
234f0b65372SJed Brown   Vec          D_loc;
235f0b65372SJed Brown   PetscScalar *d;
236f0b65372SJed Brown   PetscMemType mem_type;
237f0b65372SJed Brown 
238f0b65372SJed Brown   PetscFunctionBeginUser;
239f0b65372SJed Brown   MatShellGetContext(A, &user);
240f0b65372SJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
241f0b65372SJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
242f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2432b916ea7SJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
244f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
245f0b65372SJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
246f0b65372SJed Brown   PetscCall(VecZeroEntries(D));
247f0b65372SJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
248f0b65372SJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
249f0b65372SJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
250f0b65372SJed Brown   PetscFunctionReturn(0);
251f0b65372SJed Brown }
252f0b65372SJed Brown 
2532b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
254b107fddaSJed Brown   PetscCount ncoo;
255b107fddaSJed Brown   PetscInt  *rows, *cols;
256b107fddaSJed Brown 
257b107fddaSJed Brown   PetscFunctionBeginUser;
258b107fddaSJed Brown   if (pbdiagonal) {
259b107fddaSJed Brown     CeedSize l_size;
260b107fddaSJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
261b107fddaSJed Brown     ncoo = l_size * 5;
262b107fddaSJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
263b107fddaSJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
264b107fddaSJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
265b107fddaSJed Brown       for (PetscInt i = 0; i < 5; i++) {
266b107fddaSJed Brown         for (PetscInt j = 0; j < 5; j++) {
267b107fddaSJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
268b107fddaSJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
269b107fddaSJed Brown         }
270b107fddaSJed Brown       }
271b107fddaSJed Brown     }
272b107fddaSJed Brown   } else {
2732b916ea7SJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
274b107fddaSJed Brown   }
275b107fddaSJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
276b107fddaSJed Brown   free(rows);
277b107fddaSJed Brown   free(cols);
278b107fddaSJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
279b107fddaSJed Brown   PetscFunctionReturn(0);
280b107fddaSJed Brown }
281b107fddaSJed Brown 
2822b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
283b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
284b107fddaSJed Brown   const PetscScalar *values;
285b107fddaSJed Brown   MatType            mat_type;
286b107fddaSJed Brown 
287b107fddaSJed Brown   PetscFunctionBeginUser;
288b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
2892b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
290b107fddaSJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
2912b916ea7SJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
292b107fddaSJed Brown   } else {
293b107fddaSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
294b107fddaSJed Brown   }
295b107fddaSJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
296b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
297b107fddaSJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
298b107fddaSJed Brown   PetscFunctionReturn(0);
299b107fddaSJed Brown }
300b107fddaSJed Brown 
3012b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
302f0b65372SJed Brown   User      user = *(User *)user_data;
30304855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
304f0b65372SJed Brown   PetscFunctionBeginUser;
3052b916ea7SJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorContextSetDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
30604855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
307f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3082b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
309f0b65372SJed Brown   if (!user->matrices_set_up) {
310f0b65372SJed Brown     if (J_is_shell) {
311f0b65372SJed Brown       PetscCall(MatShellSetContext(J, user));
3122b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3132b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
314f0b65372SJed Brown       PetscCall(MatSetUp(J));
315f0b65372SJed Brown     }
316f0b65372SJed Brown     if (!J_pre_is_shell) {
3172b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
318b107fddaSJed Brown     }
31904855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
320b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
321b107fddaSJed Brown     }
322f0b65372SJed Brown     user->matrices_set_up = true;
323f0b65372SJed Brown   }
324f0b65372SJed Brown   if (!J_pre_is_shell) {
3252b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
326f0b65372SJed Brown   }
32704855949SJed Brown   if (user->coo_values_amat) {
32804855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
32904855949SJed Brown   } else if (J_is_mffd) {
33004855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
33104855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
33204855949SJed Brown   }
333f0b65372SJed Brown   PetscFunctionReturn(0);
334f0b65372SJed Brown }
335f0b65372SJed Brown 
3362b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
337a515125bSLeila Ghaffari   Vec         Q_loc;
338a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
339a515125bSLeila Ghaffari   PetscViewer viewer;
340a515125bSLeila Ghaffari   PetscFunctionBeginUser;
341a515125bSLeila Ghaffari 
342852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
343a515125bSLeila Ghaffari     // Set up output
3447538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
3457538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3467538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3477538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
348a515125bSLeila Ghaffari 
349a515125bSLeila Ghaffari     // Output
350852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3517538d537SJames Wright 
3522b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
3537538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
3547538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
355a515125bSLeila Ghaffari     if (user->dm_viz) {
356a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
357a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
358a515125bSLeila Ghaffari       PetscViewer viewer_refined;
359a515125bSLeila Ghaffari 
3607538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
3617538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
3627538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
3637538d537SJames Wright 
3647538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
3657538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3662b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
3677538d537SJames Wright 
368852e5969SJed Brown       PetscCall(
369852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3707538d537SJames Wright 
3712b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
3727538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
3737538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
3747538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
3757538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
376a515125bSLeila Ghaffari     }
3777538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
378852e5969SJed Brown   }
379a515125bSLeila Ghaffari 
380a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
38191a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
382852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
38391a36801SJames Wright   } else {
3842b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
38591a36801SJames Wright   }
3862b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
3877538d537SJames Wright 
3889293eaa1SJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
3899293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
3909293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3919293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
3929293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
3937538d537SJames Wright   PetscCall(VecView(Q, viewer));
3947538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
3957538d537SJames Wright   PetscFunctionReturn(0);
3967538d537SJames Wright }
3977538d537SJames Wright 
3987538d537SJames Wright // User provided TS Monitor
3992b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
4007538d537SJames Wright   User user = ctx;
4017538d537SJames Wright   PetscFunctionBeginUser;
4027538d537SJames Wright 
403852e5969SJed Brown   // Print every 'checkpoint_interval' steps
404c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
405c539088bSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
406c539088bSJames Wright     PetscFunctionReturn(0);
4077538d537SJames Wright 
4087538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
409a515125bSLeila Ghaffari 
410a515125bSLeila Ghaffari   PetscFunctionReturn(0);
411a515125bSLeila Ghaffari }
412a515125bSLeila Ghaffari 
413a515125bSLeila Ghaffari // TS: Create, setup, and solve
4142b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
415a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
416a515125bSLeila Ghaffari   TSAdapt     adapt;
417a515125bSLeila Ghaffari   PetscScalar final_time;
418a515125bSLeila Ghaffari   PetscFunctionBeginUser;
419a515125bSLeila Ghaffari 
4202b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4212b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
422a515125bSLeila Ghaffari   if (phys->implicit) {
4232b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
424a515125bSLeila Ghaffari     if (user->op_ifunction) {
4252b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
426a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4272b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
428a515125bSLeila Ghaffari     }
429f0b65372SJed Brown     if (user->op_ijacobian) {
4302b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
431b107fddaSJed Brown       if (app_ctx->amat_type) {
432b107fddaSJed Brown         Mat Pmat, Amat;
4332b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4342b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4352b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4362b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4372b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4382b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
439b107fddaSJed Brown       }
440f0b65372SJed Brown     }
441a515125bSLeila Ghaffari   } else {
4422b916ea7SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4432b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4442b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4452b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
446a515125bSLeila Ghaffari   }
4472b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4482b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
449f0784ed3SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4502b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4512b916ea7SJeremy L Thompson   if (app_ctx->test_mode) {
4522b916ea7SJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4532b916ea7SJeremy L Thompson   }
4542b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4552b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4562b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
457e2f84137SJeremy L Thompson   user->time = -1.0;  // require all BCs and ctx to be updated
458e2f84137SJeremy L Thompson   user->dt   = -1.0;
459a515125bSLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
460a515125bSLeila Ghaffari     if (!app_ctx->test_mode) {
4612b916ea7SJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
462a515125bSLeila Ghaffari     }
463a515125bSLeila Ghaffari   } else {  // continue from time of last output
464a515125bSLeila Ghaffari     PetscInt    count;
465a515125bSLeila Ghaffari     PetscViewer viewer;
4662b916ea7SJeremy L Thompson 
4679293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4682b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4699293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4702b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4719293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4729293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4739293eaa1SJed Brown     }
4749293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
47574a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
476a515125bSLeila Ghaffari   }
477a515125bSLeila Ghaffari   if (!app_ctx->test_mode) {
4782b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
479a515125bSLeila Ghaffari   }
480a515125bSLeila Ghaffari 
481a515125bSLeila Ghaffari   // Solve
48274a6f4ddSJed Brown   PetscReal start_time;
48374a6f4ddSJed Brown   PetscInt  start_step;
4842b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
48574a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
48691982731SJeremy L Thompson 
48791982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
48891982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
48974a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
49091982731SJeremy L Thompson   if (PetscPreLoadingOn) {
49191982731SJeremy L Thompson     // LCOV_EXCL_START
49291982731SJeremy L Thompson     SNES      snes;
49391982731SJeremy L Thompson     Vec       Q_preload;
49491982731SJeremy L Thompson     PetscReal rtol;
49591982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
49691982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
49791982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
49891982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
4992b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
50091982731SJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
50191982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5022b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
50391982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
50491982731SJeremy L Thompson     // LCOV_EXCL_STOP
50591982731SJeremy L Thompson   } else {
5062b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5072b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
50891982731SJeremy L Thompson   }
50991982731SJeremy L Thompson   PetscPreLoadEnd();
51091982731SJeremy L Thompson 
51191982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
512a515125bSLeila Ghaffari   *f_time = final_time;
51391982731SJeremy L Thompson 
514a515125bSLeila Ghaffari   if (!app_ctx->test_mode) {
515852e5969SJed Brown     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
5167538d537SJames Wright       PetscInt step_no;
5177538d537SJames Wright       PetscCall(TSGetStepNumber(*ts, &step_no));
5187538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
5197538d537SJames Wright     }
5207538d537SJames Wright 
52191982731SJeremy L Thompson     PetscLogEvent stage_id;
52291982731SJeremy L Thompson     PetscStageLog stage_log;
52391982731SJeremy L Thompson 
52491982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
52591982731SJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
5262b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
527a515125bSLeila Ghaffari   }
528a515125bSLeila Ghaffari   PetscFunctionReturn(0);
529a515125bSLeila Ghaffari }
530