xref: /libCEED/examples/fluids/src/setupts.c (revision 7b39487de4ab72ae5651458a5c845e3d06436385)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #include "../navierstokes.h"
1277841947SLeila Ghaffari #include "../qfunctions/mass.h"
1377841947SLeila Ghaffari 
1477841947SLeila Ghaffari // Compute mass matrix for explicit scheme
1577841947SLeila Ghaffari PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data,
1677841947SLeila Ghaffari                                        Vec M) {
1777841947SLeila Ghaffari   Vec            M_loc;
1877841947SLeila Ghaffari   CeedQFunction  qf_mass;
1977841947SLeila Ghaffari   CeedOperator   op_mass;
2077841947SLeila Ghaffari   CeedVector     m_ceed, ones_vec;
2177841947SLeila Ghaffari   CeedInt        num_comp_q, q_data_size;
2277841947SLeila Ghaffari   PetscErrorCode ierr;
2377841947SLeila Ghaffari   PetscFunctionBeginUser;
2477841947SLeila Ghaffari 
2577841947SLeila Ghaffari   // CEED Restriction
2677841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
2777841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
2877841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
2977841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
3077841947SLeila Ghaffari   CeedVectorSetValue(ones_vec, 1.0);
3177841947SLeila Ghaffari 
3277841947SLeila Ghaffari   // CEED QFunction
3377841947SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
3477841947SLeila Ghaffari   CeedQFunctionAddInput(qf_mass, "q", num_comp_q, CEED_EVAL_INTERP);
35a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE);
3677841947SLeila Ghaffari   CeedQFunctionAddOutput(qf_mass, "v", num_comp_q, CEED_EVAL_INTERP);
3777841947SLeila Ghaffari 
3877841947SLeila Ghaffari   // CEED Operator
3977841947SLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
4077841947SLeila Ghaffari   CeedOperatorSetField(op_mass, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
4177841947SLeila Ghaffari                        CEED_VECTOR_ACTIVE);
42a61c78d6SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i,
4377841947SLeila Ghaffari                        CEED_BASIS_COLLOCATED, ceed_data->q_data);
4477841947SLeila Ghaffari   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
4577841947SLeila Ghaffari                        CEED_VECTOR_ACTIVE);
4677841947SLeila Ghaffari 
4777841947SLeila Ghaffari   // Place PETSc vector in CEED vector
4877841947SLeila Ghaffari   CeedScalar *m;
4977841947SLeila Ghaffari   PetscMemType m_mem_type;
5077841947SLeila Ghaffari   ierr = DMGetLocalVector(dm, &M_loc); CHKERRQ(ierr);
5177841947SLeila Ghaffari   ierr = VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type);
5277841947SLeila Ghaffari   CHKERRQ(ierr);
5377841947SLeila Ghaffari   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
5477841947SLeila Ghaffari 
5577841947SLeila Ghaffari   // Apply CEED Operator
5677841947SLeila Ghaffari   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
5777841947SLeila Ghaffari 
5877841947SLeila Ghaffari   // Restore vectors
5977841947SLeila Ghaffari   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
6077841947SLeila Ghaffari   ierr = VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m);
6177841947SLeila Ghaffari   CHKERRQ(ierr);
6277841947SLeila Ghaffari 
6377841947SLeila Ghaffari   // Local-to-Global
6477841947SLeila Ghaffari   ierr = VecZeroEntries(M); CHKERRQ(ierr);
6577841947SLeila Ghaffari   ierr = DMLocalToGlobal(dm, M_loc, ADD_VALUES, M); CHKERRQ(ierr);
6677841947SLeila Ghaffari   ierr = DMRestoreLocalVector(dm, &M_loc); CHKERRQ(ierr);
6777841947SLeila Ghaffari 
6877841947SLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
6977841947SLeila Ghaffari   ierr = VecReciprocal(M); CHKERRQ(ierr);
7077841947SLeila Ghaffari 
7177841947SLeila Ghaffari   // Cleanup
7277841947SLeila Ghaffari   CeedVectorDestroy(&ones_vec);
7377841947SLeila Ghaffari   CeedVectorDestroy(&m_ceed);
7477841947SLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
7577841947SLeila Ghaffari   CeedOperatorDestroy(&op_mass);
7677841947SLeila Ghaffari 
7777841947SLeila Ghaffari   PetscFunctionReturn(0);
7877841947SLeila Ghaffari }
7977841947SLeila Ghaffari 
8077841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
8177841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
8277841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
8377841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
8477841947SLeila Ghaffari   User           user = *(User *)user_data;
8577841947SLeila Ghaffari   PetscScalar    *q, *g;
865e82a6e1SJeremy L Thompson   Vec            Q_loc = user->Q_loc, G_loc;
8777841947SLeila Ghaffari   PetscMemType   q_mem_type, g_mem_type;
8877841947SLeila Ghaffari   PetscErrorCode ierr;
8977841947SLeila Ghaffari   PetscFunctionBeginUser;
9077841947SLeila Ghaffari 
915e82a6e1SJeremy L Thompson   // Get local vector
925e82a6e1SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
935e82a6e1SJeremy L Thompson 
945e82a6e1SJeremy L Thompson   // Update time dependent data
955e82a6e1SJeremy L Thompson   if (user->time != t) {
965e82a6e1SJeremy L Thompson     ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t,
975e82a6e1SJeremy L Thompson                                       NULL, NULL, NULL); CHKERRQ(ierr);
985e82a6e1SJeremy L Thompson     if (user->phys->solution_time_label) {
9911436a05SJames Wright       CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t);
1005e82a6e1SJeremy L Thompson     }
1015e82a6e1SJeremy L Thompson     user->time = t;
1025e82a6e1SJeremy L Thompson   }
10388626eedSJames Wright   if (user->phys->timestep_size_label) {
10488626eedSJames Wright     PetscScalar dt;
10588626eedSJames Wright     ierr = TSGetTimeStep(ts, &dt); CHKERRQ(ierr);
1065e82a6e1SJeremy L Thompson     if (user->dt != dt) {
10788626eedSJames Wright       CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label,
10888626eedSJames Wright                                    &dt);
1095e82a6e1SJeremy L Thompson       user->dt = dt;
1105e82a6e1SJeremy L Thompson     }
11188626eedSJames Wright   }
11277841947SLeila Ghaffari 
11377841947SLeila Ghaffari   // Global-to-local
11477841947SLeila Ghaffari   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
11577841947SLeila Ghaffari 
11677841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
11777841947SLeila Ghaffari   ierr = VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type);
11877841947SLeila Ghaffari   CHKERRQ(ierr);
11977841947SLeila Ghaffari   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
12077841947SLeila Ghaffari   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
12177841947SLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
12277841947SLeila Ghaffari 
12377841947SLeila Ghaffari   // Apply CEED operator
12477841947SLeila Ghaffari   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed,
12577841947SLeila Ghaffari                     CEED_REQUEST_IMMEDIATE);
12677841947SLeila Ghaffari 
12777841947SLeila Ghaffari   // Restore vectors
12877841947SLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
12977841947SLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
13077841947SLeila Ghaffari   ierr = VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q);
13177841947SLeila Ghaffari   CHKERRQ(ierr);
13277841947SLeila Ghaffari   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
13377841947SLeila Ghaffari 
13477841947SLeila Ghaffari   // Local-to-Global
13577841947SLeila Ghaffari   ierr = VecZeroEntries(G); CHKERRQ(ierr);
13677841947SLeila Ghaffari   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
13777841947SLeila Ghaffari 
13877841947SLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
13977841947SLeila Ghaffari   ierr = VecPointwiseMult(G, G, user->M); CHKERRQ(ierr);
14077841947SLeila Ghaffari 
14177841947SLeila Ghaffari   // Restore vectors
14277841947SLeila Ghaffari   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
14377841947SLeila Ghaffari 
14477841947SLeila Ghaffari   PetscFunctionReturn(0);
14577841947SLeila Ghaffari }
14677841947SLeila Ghaffari 
14777841947SLeila Ghaffari // Implicit time-stepper function setup
14877841947SLeila Ghaffari PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G,
14977841947SLeila Ghaffari                             void *user_data) {
15077841947SLeila Ghaffari   User              user = *(User *)user_data;
15177841947SLeila Ghaffari   const PetscScalar *q, *q_dot;
15277841947SLeila Ghaffari   PetscScalar       *g;
1535e82a6e1SJeremy L Thompson   Vec               Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
15477841947SLeila Ghaffari   PetscMemType      q_mem_type, q_dot_mem_type, g_mem_type;
15577841947SLeila Ghaffari   PetscErrorCode    ierr;
15677841947SLeila Ghaffari   PetscFunctionBeginUser;
15777841947SLeila Ghaffari 
1585e82a6e1SJeremy L Thompson   // Get local vectors
1595e82a6e1SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
1605e82a6e1SJeremy L Thompson 
1615e82a6e1SJeremy L Thompson   // Update time dependent data
1625e82a6e1SJeremy L Thompson   if (user->time != t) {
1635e82a6e1SJeremy L Thompson     ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t,
1645e82a6e1SJeremy L Thompson                                       NULL, NULL, NULL); CHKERRQ(ierr);
1655e82a6e1SJeremy L Thompson     if (user->phys->solution_time_label) {
16611436a05SJames Wright       CeedOperatorContextSetDouble(user->op_ifunction,
16711436a05SJames Wright                                    user->phys->solution_time_label, &t);
1685e82a6e1SJeremy L Thompson     }
1695e82a6e1SJeremy L Thompson     user->time = t;
1705e82a6e1SJeremy L Thompson   }
17188626eedSJames Wright   if (user->phys->timestep_size_label) {
17288626eedSJames Wright     PetscScalar dt;
17388626eedSJames Wright     ierr = TSGetTimeStep(ts, &dt); CHKERRQ(ierr);
1745e82a6e1SJeremy L Thompson     if (user->dt != dt) {
17588626eedSJames Wright       CeedOperatorContextSetDouble(user->op_ifunction,
17688626eedSJames Wright                                    user->phys->timestep_size_label, &dt);
1775e82a6e1SJeremy L Thompson       user->dt = dt;
1785e82a6e1SJeremy L Thompson     }
17988626eedSJames Wright   }
18077841947SLeila Ghaffari 
18177841947SLeila Ghaffari   // Global-to-local
18277841947SLeila Ghaffari   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
18377841947SLeila Ghaffari   ierr = DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc);
18477841947SLeila Ghaffari   CHKERRQ(ierr);
18577841947SLeila Ghaffari 
18677841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
18777841947SLeila Ghaffari   ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr);
18877841947SLeila Ghaffari   ierr = VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type);
18977841947SLeila Ghaffari   CHKERRQ(ierr);
19077841947SLeila Ghaffari   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
19177841947SLeila Ghaffari   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER,
19277841947SLeila Ghaffari                      (PetscScalar *)q);
19377841947SLeila Ghaffari   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type),
194e6225c47SLeila Ghaffari                      CEED_USE_POINTER, (PetscScalar *)q_dot);
19577841947SLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
19677841947SLeila Ghaffari 
19777841947SLeila Ghaffari   // Apply CEED operator
19877841947SLeila Ghaffari   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed,
19977841947SLeila Ghaffari                     CEED_REQUEST_IMMEDIATE);
20077841947SLeila Ghaffari 
20177841947SLeila Ghaffari   // Restore vectors
20277841947SLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
20377841947SLeila Ghaffari   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
20477841947SLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
20577841947SLeila Ghaffari   ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr);
20677841947SLeila Ghaffari   ierr = VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot); CHKERRQ(ierr);
20777841947SLeila Ghaffari   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
20877841947SLeila Ghaffari 
20977841947SLeila Ghaffari   // Local-to-Global
21077841947SLeila Ghaffari   ierr = VecZeroEntries(G); CHKERRQ(ierr);
21177841947SLeila Ghaffari   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
21277841947SLeila Ghaffari 
21377841947SLeila Ghaffari   // Restore vectors
21477841947SLeila Ghaffari   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
21577841947SLeila Ghaffari 
21677841947SLeila Ghaffari   PetscFunctionReturn(0);
21777841947SLeila Ghaffari }
21877841947SLeila Ghaffari 
219e334ad8fSJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
220e334ad8fSJed Brown   User user;
221e334ad8fSJed Brown   const PetscScalar *q;
222e334ad8fSJed Brown   PetscScalar       *g;
223e334ad8fSJed Brown   PetscMemType      q_mem_type, g_mem_type;
224e334ad8fSJed Brown   PetscErrorCode    ierr;
225e334ad8fSJed Brown   PetscFunctionBeginUser;
2265e82a6e1SJeremy L Thompson   ierr = MatShellGetContext(J, &user); CHKERRQ(ierr);
2275e82a6e1SJeremy L Thompson   Vec               Q_loc = user->Q_dot_loc, // Note - Q_dot_loc has zero BCs
2285e82a6e1SJeremy L Thompson                     G_loc;
2295e82a6e1SJeremy L Thompson 
230e334ad8fSJed Brown   // Get local vectors
231e334ad8fSJed Brown   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
232e334ad8fSJed Brown 
233e334ad8fSJed Brown   // Global-to-local
234e334ad8fSJed Brown   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
235e334ad8fSJed Brown 
236e334ad8fSJed Brown   // Place PETSc vectors in CEED vectors
237e334ad8fSJed Brown   ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr);
238e334ad8fSJed Brown   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
239e334ad8fSJed Brown   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER,
240e334ad8fSJed Brown                      (PetscScalar *)q);
241e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
242e334ad8fSJed Brown 
243e334ad8fSJed Brown   // Apply CEED operator
244e334ad8fSJed Brown   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed,
245e334ad8fSJed Brown                     CEED_REQUEST_IMMEDIATE);
246e334ad8fSJed Brown 
247e334ad8fSJed Brown   // Restore vectors
248e334ad8fSJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
249e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
250e334ad8fSJed Brown   ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr);
251e334ad8fSJed Brown   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
252e334ad8fSJed Brown 
253e334ad8fSJed Brown   // Local-to-Global
254e334ad8fSJed Brown   ierr = VecZeroEntries(G); CHKERRQ(ierr);
255e334ad8fSJed Brown   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
256e334ad8fSJed Brown 
257e334ad8fSJed Brown   // Restore vectors
258e334ad8fSJed Brown   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
259e334ad8fSJed Brown   PetscFunctionReturn(0);
260e334ad8fSJed Brown }
261e334ad8fSJed Brown 
262e334ad8fSJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
263e334ad8fSJed Brown   User user;
264e334ad8fSJed Brown   Vec D_loc;
265e334ad8fSJed Brown   PetscScalar *d;
266e334ad8fSJed Brown   PetscMemType mem_type;
267e334ad8fSJed Brown 
268e334ad8fSJed Brown   PetscFunctionBeginUser;
269e334ad8fSJed Brown   MatShellGetContext(A, &user);
270e334ad8fSJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
271e334ad8fSJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
272e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
273e334ad8fSJed Brown   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed,
274e334ad8fSJed Brown                                      CEED_REQUEST_IMMEDIATE);
275e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
276e334ad8fSJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
277e334ad8fSJed Brown   PetscCall(VecZeroEntries(D));
278e334ad8fSJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
279e334ad8fSJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
280e334ad8fSJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
281e334ad8fSJed Brown   PetscFunctionReturn(0);
282e334ad8fSJed Brown }
283e334ad8fSJed Brown 
284e334ad8fSJed Brown PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot,
285e334ad8fSJed Brown                                 PetscReal shift, Mat J, Mat J_pre,
286e334ad8fSJed Brown                                 void *user_data) {
287e334ad8fSJed Brown   User user = *(User *)user_data;
288e334ad8fSJed Brown   PetscBool J_is_shell, J_pre_is_shell;
289e334ad8fSJed Brown   PetscFunctionBeginUser;
290e334ad8fSJed Brown   if (user->phys->ijacobian_time_shift_label)
291e334ad8fSJed Brown     CeedOperatorContextSetDouble(user->op_ijacobian,
292e334ad8fSJed Brown                                  user->phys->ijacobian_time_shift_label, &shift);
293e334ad8fSJed Brown   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
294e334ad8fSJed Brown   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
295e334ad8fSJed Brown   Vec coo_vec = NULL;
296e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
297e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL,
298e334ad8fSJed Brown                                    &J_pre_is_shell));
299e334ad8fSJed Brown   if (!user->matrices_set_up) {
300e334ad8fSJed Brown     if (J_is_shell) {
301e334ad8fSJed Brown       PetscCall(MatShellSetContext(J, user));
302e334ad8fSJed Brown       PetscCall(MatShellSetOperation(J, MATOP_MULT,
303e334ad8fSJed Brown                                      (void (*)(void))MatMult_NS_IJacobian));
304e334ad8fSJed Brown       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL,
305e334ad8fSJed Brown                                      (void (*)(void))MatGetDiagonal_NS_IJacobian));
306e334ad8fSJed Brown       PetscCall(MatSetUp(J));
307e334ad8fSJed Brown     }
308e334ad8fSJed Brown     if (!J_pre_is_shell) {
309e334ad8fSJed Brown       PetscCount ncoo;
310e334ad8fSJed Brown       PetscInt *rows, *cols;
311e334ad8fSJed Brown       PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows,
312e334ad8fSJed Brown                 &cols));
313e334ad8fSJed Brown       PetscCall(MatSetPreallocationCOOLocal(J_pre, ncoo, rows, cols));
314e334ad8fSJed Brown       free(rows);
315e334ad8fSJed Brown       free(cols);
316e334ad8fSJed Brown       CeedVectorCreate(user->ceed, ncoo, &user->coo_values);
317e334ad8fSJed Brown       user->matrices_set_up = true;
318e334ad8fSJed Brown       VecCreateSeq(PETSC_COMM_WORLD, ncoo, &coo_vec);
319e334ad8fSJed Brown     }
320e334ad8fSJed Brown   }
321e334ad8fSJed Brown   if (!J_pre_is_shell) {
322e334ad8fSJed Brown     CeedMemType mem_type = CEED_MEM_HOST;
323e334ad8fSJed Brown     const PetscScalar *values;
324e334ad8fSJed Brown     MatType mat_type;
325e334ad8fSJed Brown     PetscCall(MatGetType(J_pre, &mat_type));
3269a9f72ebSJeremy L Thompson     if (strstr(mat_type, "kokkos")
3279a9f72ebSJeremy L Thompson         || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
328e334ad8fSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, user->coo_values);
329e334ad8fSJed Brown     CeedVectorGetArrayRead(user->coo_values, mem_type, &values);
330e334ad8fSJed Brown     if (coo_vec) {
331e334ad8fSJed Brown       VecPlaceArray(coo_vec, values);
332e334ad8fSJed Brown       VecViewFromOptions(coo_vec, NULL, "-coo_vec_view");
333e334ad8fSJed Brown       VecDestroy(&coo_vec);
334e334ad8fSJed Brown     }
335e334ad8fSJed Brown     PetscCall(MatSetValuesCOO(J_pre, values, INSERT_VALUES));
336e334ad8fSJed Brown     CeedVectorRestoreArrayRead(user->coo_values, &values);
337e334ad8fSJed Brown   }
338e334ad8fSJed Brown   PetscFunctionReturn(0);
339e334ad8fSJed Brown }
340e334ad8fSJed Brown 
34177841947SLeila Ghaffari // User provided TS Monitor
34277841947SLeila Ghaffari PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time,
34377841947SLeila Ghaffari                             Vec Q, void *ctx) {
34477841947SLeila Ghaffari   User           user = ctx;
34577841947SLeila Ghaffari   Vec            Q_loc;
34677841947SLeila Ghaffari   char           file_path[PETSC_MAX_PATH_LEN];
34777841947SLeila Ghaffari   PetscViewer    viewer;
34877841947SLeila Ghaffari   PetscErrorCode ierr;
34977841947SLeila Ghaffari   PetscFunctionBeginUser;
35077841947SLeila Ghaffari 
35177841947SLeila Ghaffari   // Print every 'output_freq' steps
35277841947SLeila Ghaffari   if (step_no % user->app_ctx->output_freq != 0)
35377841947SLeila Ghaffari     PetscFunctionReturn(0);
35477841947SLeila Ghaffari 
35577841947SLeila Ghaffari   // Set up output
35677841947SLeila Ghaffari   ierr = DMGetLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
35777841947SLeila Ghaffari   ierr = PetscObjectSetName((PetscObject)Q_loc, "StateVec"); CHKERRQ(ierr);
35877841947SLeila Ghaffari   ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr);
35977841947SLeila Ghaffari   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
36077841947SLeila Ghaffari 
36177841947SLeila Ghaffari   // Output
36208140895SJed Brown   ierr = PetscSNPrintf(file_path, sizeof file_path,
36308140895SJed Brown                        "%s/ns-%03" PetscInt_FMT ".vtu",
36477841947SLeila Ghaffari                        user->app_ctx->output_dir, step_no + user->app_ctx->cont_steps);
36577841947SLeila Ghaffari   CHKERRQ(ierr);
36677841947SLeila Ghaffari   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path,
36777841947SLeila Ghaffari                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
36877841947SLeila Ghaffari   ierr = VecView(Q_loc, viewer); CHKERRQ(ierr);
36977841947SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
37077841947SLeila Ghaffari   if (user->dm_viz) {
37177841947SLeila Ghaffari     Vec         Q_refined, Q_refined_loc;
37277841947SLeila Ghaffari     char        file_path_refined[PETSC_MAX_PATH_LEN];
37377841947SLeila Ghaffari     PetscViewer viewer_refined;
37477841947SLeila Ghaffari 
37577841947SLeila Ghaffari     ierr = DMGetGlobalVector(user->dm_viz, &Q_refined); CHKERRQ(ierr);
37677841947SLeila Ghaffari     ierr = DMGetLocalVector(user->dm_viz, &Q_refined_loc); CHKERRQ(ierr);
37777841947SLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)Q_refined_loc, "Refined");
37877841947SLeila Ghaffari     CHKERRQ(ierr);
37977841947SLeila Ghaffari     ierr = MatInterpolate(user->interp_viz, Q, Q_refined); CHKERRQ(ierr);
38077841947SLeila Ghaffari     ierr = VecZeroEntries(Q_refined_loc); CHKERRQ(ierr);
38177841947SLeila Ghaffari     ierr = DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc);
38277841947SLeila Ghaffari     CHKERRQ(ierr);
38377841947SLeila Ghaffari     ierr = PetscSNPrintf(file_path_refined, sizeof file_path_refined,
38408140895SJed Brown                          "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir,
385e6225c47SLeila Ghaffari                          step_no + user->app_ctx->cont_steps);
38677841947SLeila Ghaffari     CHKERRQ(ierr);
38777841947SLeila Ghaffari     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined),
388e6225c47SLeila Ghaffari                               file_path_refined, FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
38977841947SLeila Ghaffari     ierr = VecView(Q_refined_loc, viewer_refined); CHKERRQ(ierr);
39077841947SLeila Ghaffari     ierr = DMRestoreLocalVector(user->dm_viz, &Q_refined_loc); CHKERRQ(ierr);
39177841947SLeila Ghaffari     ierr = DMRestoreGlobalVector(user->dm_viz, &Q_refined); CHKERRQ(ierr);
39277841947SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
39377841947SLeila Ghaffari   }
39477841947SLeila Ghaffari   ierr = DMRestoreLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
39577841947SLeila Ghaffari 
39677841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
39777841947SLeila Ghaffari   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin",
39877841947SLeila Ghaffari                        user->app_ctx->output_dir); CHKERRQ(ierr);
39977841947SLeila Ghaffari   ierr = PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer);
40077841947SLeila Ghaffari   CHKERRQ(ierr);
40177841947SLeila Ghaffari   ierr = VecView(Q, viewer); CHKERRQ(ierr);
40277841947SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
40377841947SLeila Ghaffari 
40477841947SLeila Ghaffari   // Save time stamp
40577841947SLeila Ghaffari   // Dimensionalize time back
40677841947SLeila Ghaffari   time /= user->units->second;
40777841947SLeila Ghaffari   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin",
40877841947SLeila Ghaffari                        user->app_ctx->output_dir); CHKERRQ(ierr);
40977841947SLeila Ghaffari   ierr = PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer);
41077841947SLeila Ghaffari   CHKERRQ(ierr);
41177841947SLeila Ghaffari   #if PETSC_VERSION_GE(3,13,0)
41277841947SLeila Ghaffari   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
41377841947SLeila Ghaffari   #else
41477841947SLeila Ghaffari   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
41577841947SLeila Ghaffari   #endif
41677841947SLeila Ghaffari   CHKERRQ(ierr);
41777841947SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
41877841947SLeila Ghaffari 
41977841947SLeila Ghaffari   PetscFunctionReturn(0);
42077841947SLeila Ghaffari }
42177841947SLeila Ghaffari 
42277841947SLeila Ghaffari // TS: Create, setup, and solve
42377841947SLeila Ghaffari PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys,
42477841947SLeila Ghaffari                           Vec *Q, PetscScalar *f_time, TS *ts) {
42577841947SLeila Ghaffari   MPI_Comm       comm = user->comm;
42677841947SLeila Ghaffari   TSAdapt        adapt;
42777841947SLeila Ghaffari   PetscScalar    final_time;
42877841947SLeila Ghaffari   PetscErrorCode ierr;
42977841947SLeila Ghaffari   PetscFunctionBeginUser;
43077841947SLeila Ghaffari 
43177841947SLeila Ghaffari   ierr = TSCreate(comm, ts); CHKERRQ(ierr);
43277841947SLeila Ghaffari   ierr = TSSetDM(*ts, dm); CHKERRQ(ierr);
43377841947SLeila Ghaffari   if (phys->implicit) {
43477841947SLeila Ghaffari     ierr = TSSetType(*ts, TSBDF); CHKERRQ(ierr);
43577841947SLeila Ghaffari     if (user->op_ifunction) {
43677841947SLeila Ghaffari       ierr = TSSetIFunction(*ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
43777841947SLeila Ghaffari     } else { // Implicit integrators can fall back to using an RHSFunction
43877841947SLeila Ghaffari       ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
43977841947SLeila Ghaffari     }
440e334ad8fSJed Brown     if (user->op_ijacobian) {
441e334ad8fSJed Brown       ierr = DMTSSetIJacobian(dm, FormIJacobian_NS, &user); CHKERRQ(ierr);
442e334ad8fSJed Brown     }
44377841947SLeila Ghaffari   } else {
44477841947SLeila Ghaffari     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
44577841947SLeila Ghaffari                                  "Problem does not provide RHSFunction");
44677841947SLeila Ghaffari     ierr = TSSetType(*ts, TSRK); CHKERRQ(ierr);
44777841947SLeila Ghaffari     ierr = TSRKSetType(*ts, TSRK5F); CHKERRQ(ierr);
44877841947SLeila Ghaffari     ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
44977841947SLeila Ghaffari   }
45077841947SLeila Ghaffari   ierr = TSSetMaxTime(*ts, 500. * user->units->second); CHKERRQ(ierr);
45177841947SLeila Ghaffari   ierr = TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
45277841947SLeila Ghaffari   ierr = TSSetTimeStep(*ts, 1.e-2 * user->units->second); CHKERRQ(ierr);
45377841947SLeila Ghaffari   if (app_ctx->test_mode) {ierr = TSSetMaxSteps(*ts, 10); CHKERRQ(ierr);}
45477841947SLeila Ghaffari   ierr = TSGetAdapt(*ts, &adapt); CHKERRQ(ierr);
45577841947SLeila Ghaffari   ierr = TSAdaptSetStepLimits(adapt,
45677841947SLeila Ghaffari                               1.e-12 * user->units->second,
45777841947SLeila Ghaffari                               1.e2 * user->units->second); CHKERRQ(ierr);
45877841947SLeila Ghaffari   ierr = TSSetFromOptions(*ts); CHKERRQ(ierr);
4595e82a6e1SJeremy L Thompson   user->time = -1.0; // require all BCs and ctx to be updated
4605e82a6e1SJeremy L Thompson   user->dt   = -1.0;
46177841947SLeila Ghaffari   if (!app_ctx->cont_steps) { // print initial condition
46277841947SLeila Ghaffari     if (!app_ctx->test_mode) {
46377841947SLeila Ghaffari       ierr = TSMonitor_NS(*ts, 0, 0., *Q, user); CHKERRQ(ierr);
46477841947SLeila Ghaffari     }
46577841947SLeila Ghaffari   } else { // continue from time of last output
46677841947SLeila Ghaffari     PetscReal   time;
46777841947SLeila Ghaffari     PetscInt    count;
46877841947SLeila Ghaffari     PetscViewer viewer;
46977841947SLeila Ghaffari     char        file_path[PETSC_MAX_PATH_LEN];
47077841947SLeila Ghaffari     ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin",
47177841947SLeila Ghaffari                          app_ctx->output_dir); CHKERRQ(ierr);
47277841947SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, file_path, FILE_MODE_READ, &viewer);
47377841947SLeila Ghaffari     CHKERRQ(ierr);
47477841947SLeila Ghaffari     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
47577841947SLeila Ghaffari     CHKERRQ(ierr);
47677841947SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
47777841947SLeila Ghaffari     ierr = TSSetTime(*ts, time * user->units->second); CHKERRQ(ierr);
47877841947SLeila Ghaffari   }
47977841947SLeila Ghaffari   if (!app_ctx->test_mode) {
48077841947SLeila Ghaffari     ierr = TSMonitorSet(*ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
48177841947SLeila Ghaffari   }
48277841947SLeila Ghaffari 
48377841947SLeila Ghaffari   // Solve
484*7b39487dSJeremy L Thompson   PetscScalar start_time;
485*7b39487dSJeremy L Thompson   ierr = TSGetTime(*ts, &start_time); CHKERRQ(ierr);
486*7b39487dSJeremy L Thompson 
487*7b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
488*7b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
489*7b39487dSJeremy L Thompson   PetscCall(TSSetStepNumber(*ts, 0));
490*7b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
491*7b39487dSJeremy L Thompson     // LCOV_EXCL_START
492*7b39487dSJeremy L Thompson     SNES      snes;
493*7b39487dSJeremy L Thompson     Vec       Q_preload;
494*7b39487dSJeremy L Thompson     PetscReal rtol;
495*7b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
496*7b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
497*7b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
498*7b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
499*7b39487dSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT,
500*7b39487dSJeremy L Thompson                                 PETSC_DEFAULT, PETSC_DEFAULT));
501*7b39487dSJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
502*7b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
503*7b39487dSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT,
504*7b39487dSJeremy L Thompson                                 PETSC_DEFAULT, PETSC_DEFAULT));
505*7b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
506*7b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
507*7b39487dSJeremy L Thompson   } else {
50877841947SLeila Ghaffari     ierr = PetscBarrier((PetscObject) *ts); CHKERRQ(ierr);
50977841947SLeila Ghaffari     ierr = TSSolve(*ts, *Q); CHKERRQ(ierr);
510*7b39487dSJeremy L Thompson   }
511*7b39487dSJeremy L Thompson   PetscPreLoadEnd();
512*7b39487dSJeremy L Thompson 
513*7b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
51477841947SLeila Ghaffari   *f_time = final_time;
515*7b39487dSJeremy L Thompson 
51677841947SLeila Ghaffari   if (!app_ctx->test_mode) {
517*7b39487dSJeremy L Thompson     PetscLogEvent stage_id;
518*7b39487dSJeremy L Thompson     PetscStageLog stage_log;
519*7b39487dSJeremy L Thompson 
520*7b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
521*7b39487dSJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
522*7b39487dSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD,
52377841947SLeila Ghaffari                           "Time taken for solution (sec): %g\n",
524*7b39487dSJeremy L Thompson                           stage_log->stageInfo[stage_id].perfInfo.time));
52577841947SLeila Ghaffari   }
52677841947SLeila Ghaffari   PetscFunctionReturn(0);
52777841947SLeila Ghaffari }
528