xref: /libCEED/examples/fluids/src/setupts.c (revision ef080ff9ce83a1650979d1b767b88f0d6e3ee6ca)
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 
1377841947SLeila Ghaffari // Compute mass matrix for explicit scheme
142b730f8bSJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
1577841947SLeila Ghaffari   Vec           M_loc;
1677841947SLeila Ghaffari   CeedQFunction qf_mass;
1777841947SLeila Ghaffari   CeedOperator  op_mass;
1877841947SLeila Ghaffari   CeedVector    m_ceed, ones_vec;
1977841947SLeila Ghaffari   CeedInt       num_comp_q, q_data_size;
2077841947SLeila Ghaffari   PetscFunctionBeginUser;
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari   // CEED Restriction
2377841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
2477841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
2577841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
2677841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
2777841947SLeila Ghaffari   CeedVectorSetValue(ones_vec, 1.0);
2877841947SLeila Ghaffari 
2977841947SLeila Ghaffari   // CEED QFunction
30*ef080ff9SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
3177841947SLeila Ghaffari 
3277841947SLeila Ghaffari   // CEED Operator
3377841947SLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
34*ef080ff9SJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
352b730f8bSJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
362b730f8bSJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
3777841947SLeila Ghaffari 
3877841947SLeila Ghaffari   // Place PETSc vector in CEED vector
3977841947SLeila Ghaffari   CeedScalar  *m;
4077841947SLeila Ghaffari   PetscMemType m_mem_type;
412b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &M_loc));
422b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type));
4377841947SLeila Ghaffari   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
4477841947SLeila Ghaffari 
4577841947SLeila Ghaffari   // Apply CEED Operator
4677841947SLeila Ghaffari   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
4777841947SLeila Ghaffari 
4877841947SLeila Ghaffari   // Restore vectors
4977841947SLeila Ghaffari   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
502b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m));
5177841947SLeila Ghaffari 
5277841947SLeila Ghaffari   // Local-to-Global
532b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(M));
542b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M));
552b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &M_loc));
5677841947SLeila Ghaffari 
5777841947SLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
582b730f8bSJeremy L Thompson   PetscCall(VecReciprocal(M));
5977841947SLeila Ghaffari 
6077841947SLeila Ghaffari   // Cleanup
6177841947SLeila Ghaffari   CeedVectorDestroy(&ones_vec);
6277841947SLeila Ghaffari   CeedVectorDestroy(&m_ceed);
6377841947SLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
6477841947SLeila Ghaffari   CeedOperatorDestroy(&op_mass);
6577841947SLeila Ghaffari 
6677841947SLeila Ghaffari   PetscFunctionReturn(0);
6777841947SLeila Ghaffari }
6877841947SLeila Ghaffari 
6977841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
7077841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
7177841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
7277841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
7377841947SLeila Ghaffari   User         user = *(User *)user_data;
7477841947SLeila Ghaffari   PetscScalar *q, *g;
755e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
7677841947SLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
7777841947SLeila Ghaffari   PetscFunctionBeginUser;
7877841947SLeila Ghaffari 
795e82a6e1SJeremy L Thompson   // Get local vector
802b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
815e82a6e1SJeremy L Thompson 
825e82a6e1SJeremy L Thompson   // Update time dependent data
835e82a6e1SJeremy L Thompson   if (user->time != t) {
842b730f8bSJeremy L Thompson     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
855e82a6e1SJeremy L Thompson     if (user->phys->solution_time_label) {
8611436a05SJames Wright       CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t);
875e82a6e1SJeremy L Thompson     }
885e82a6e1SJeremy L Thompson     user->time = t;
895e82a6e1SJeremy L Thompson   }
9088626eedSJames Wright   if (user->phys->timestep_size_label) {
9188626eedSJames Wright     PetscScalar dt;
922b730f8bSJeremy L Thompson     PetscCall(TSGetTimeStep(ts, &dt));
935e82a6e1SJeremy L Thompson     if (user->dt != dt) {
942b730f8bSJeremy L Thompson       CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label, &dt);
955e82a6e1SJeremy L Thompson       user->dt = dt;
965e82a6e1SJeremy L Thompson     }
9788626eedSJames Wright   }
9877841947SLeila Ghaffari 
9977841947SLeila Ghaffari   // Global-to-local
1002b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
10177841947SLeila Ghaffari 
10277841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
1032b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type));
1042b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
10577841947SLeila Ghaffari   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
10677841947SLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
10777841947SLeila Ghaffari 
10877841947SLeila Ghaffari   // Apply CEED operator
1092b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
11077841947SLeila Ghaffari 
11177841947SLeila Ghaffari   // Restore vectors
11277841947SLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
11377841947SLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
1142b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q));
1152b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
11677841947SLeila Ghaffari 
11777841947SLeila Ghaffari   // Local-to-Global
1182b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
1192b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
12077841947SLeila Ghaffari 
12177841947SLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1222b730f8bSJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
12377841947SLeila Ghaffari 
12477841947SLeila Ghaffari   // Restore vectors
1252b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
12677841947SLeila Ghaffari 
12777841947SLeila Ghaffari   PetscFunctionReturn(0);
12877841947SLeila Ghaffari }
12977841947SLeila Ghaffari 
13077841947SLeila Ghaffari // Implicit time-stepper function setup
1312b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
13277841947SLeila Ghaffari   User               user = *(User *)user_data;
13377841947SLeila Ghaffari   const PetscScalar *q, *q_dot;
13477841947SLeila Ghaffari   PetscScalar       *g;
1355e82a6e1SJeremy L Thompson   Vec                Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
13677841947SLeila Ghaffari   PetscMemType       q_mem_type, q_dot_mem_type, g_mem_type;
13777841947SLeila Ghaffari   PetscFunctionBeginUser;
13877841947SLeila Ghaffari 
1395e82a6e1SJeremy L Thompson   // Get local vectors
1402b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
1415e82a6e1SJeremy L Thompson 
1425e82a6e1SJeremy L Thompson   // Update time dependent data
1435e82a6e1SJeremy L Thompson   if (user->time != t) {
1442b730f8bSJeremy L Thompson     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
1455e82a6e1SJeremy L Thompson     if (user->phys->solution_time_label) {
1462b730f8bSJeremy L Thompson       CeedOperatorContextSetDouble(user->op_ifunction, user->phys->solution_time_label, &t);
1475e82a6e1SJeremy L Thompson     }
1485e82a6e1SJeremy L Thompson     user->time = t;
1495e82a6e1SJeremy L Thompson   }
15088626eedSJames Wright   if (user->phys->timestep_size_label) {
15188626eedSJames Wright     PetscScalar dt;
1522b730f8bSJeremy L Thompson     PetscCall(TSGetTimeStep(ts, &dt));
1535e82a6e1SJeremy L Thompson     if (user->dt != dt) {
1542b730f8bSJeremy L Thompson       CeedOperatorContextSetDouble(user->op_ifunction, user->phys->timestep_size_label, &dt);
1555e82a6e1SJeremy L Thompson       user->dt = dt;
1565e82a6e1SJeremy L Thompson     }
15788626eedSJames Wright   }
15877841947SLeila Ghaffari 
15977841947SLeila Ghaffari   // Global-to-local
1602b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
1612b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
16277841947SLeila Ghaffari 
16377841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
1642b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
1652b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type));
1662b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
1672b730f8bSJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
1682b730f8bSJeremy L Thompson   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), CEED_USE_POINTER, (PetscScalar *)q_dot);
16977841947SLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
17077841947SLeila Ghaffari 
17177841947SLeila Ghaffari   // Apply CEED operator
1722b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
17377841947SLeila Ghaffari 
17477841947SLeila Ghaffari   // Restore vectors
17577841947SLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
17677841947SLeila Ghaffari   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
17777841947SLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
1782b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
1792b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot));
1802b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
18177841947SLeila Ghaffari 
18277841947SLeila Ghaffari   // Local-to-Global
1832b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
1842b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
18577841947SLeila Ghaffari 
18677841947SLeila Ghaffari   // Restore vectors
1872b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
18877841947SLeila Ghaffari 
18977841947SLeila Ghaffari   PetscFunctionReturn(0);
19077841947SLeila Ghaffari }
19177841947SLeila Ghaffari 
192e334ad8fSJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
193e334ad8fSJed Brown   User               user;
194e334ad8fSJed Brown   const PetscScalar *q;
195e334ad8fSJed Brown   PetscScalar       *g;
196e334ad8fSJed Brown   PetscMemType       q_mem_type, g_mem_type;
197e334ad8fSJed Brown   PetscFunctionBeginUser;
1982b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
1995e82a6e1SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
2005e82a6e1SJeremy L Thompson       G_loc;
2015e82a6e1SJeremy L Thompson 
202e334ad8fSJed Brown   // Get local vectors
2032b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
204e334ad8fSJed Brown 
205e334ad8fSJed Brown   // Global-to-local
2062b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
207e334ad8fSJed Brown 
208e334ad8fSJed Brown   // Place PETSc vectors in CEED vectors
2092b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2102b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2112b730f8bSJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
212e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
213e334ad8fSJed Brown 
214e334ad8fSJed Brown   // Apply CEED operator
2152b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
216e334ad8fSJed Brown 
217e334ad8fSJed Brown   // Restore vectors
218e334ad8fSJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
219e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2202b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2212b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
222e334ad8fSJed Brown 
223e334ad8fSJed Brown   // Local-to-Global
2242b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2252b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
226e334ad8fSJed Brown 
227e334ad8fSJed Brown   // Restore vectors
2282b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
229e334ad8fSJed Brown   PetscFunctionReturn(0);
230e334ad8fSJed Brown }
231e334ad8fSJed Brown 
232e334ad8fSJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
233e334ad8fSJed Brown   User         user;
234e334ad8fSJed Brown   Vec          D_loc;
235e334ad8fSJed Brown   PetscScalar *d;
236e334ad8fSJed Brown   PetscMemType mem_type;
237e334ad8fSJed Brown 
238e334ad8fSJed Brown   PetscFunctionBeginUser;
239e334ad8fSJed Brown   MatShellGetContext(A, &user);
240e334ad8fSJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
241e334ad8fSJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
242e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2432b730f8bSJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
244e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
245e334ad8fSJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
246e334ad8fSJed Brown   PetscCall(VecZeroEntries(D));
247e334ad8fSJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
248e334ad8fSJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
249e334ad8fSJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
250e334ad8fSJed Brown   PetscFunctionReturn(0);
251e334ad8fSJed Brown }
252e334ad8fSJed Brown 
2532b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
254544be873SJed Brown   PetscCount ncoo;
255544be873SJed Brown   PetscInt  *rows, *cols;
256544be873SJed Brown 
257544be873SJed Brown   PetscFunctionBeginUser;
258544be873SJed Brown   if (pbdiagonal) {
259544be873SJed Brown     CeedSize l_size;
260544be873SJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
261544be873SJed Brown     ncoo = l_size * 5;
262544be873SJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
263544be873SJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
264544be873SJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
265544be873SJed Brown       for (PetscInt i = 0; i < 5; i++) {
266544be873SJed Brown         for (PetscInt j = 0; j < 5; j++) {
267544be873SJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
268544be873SJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
269544be873SJed Brown         }
270544be873SJed Brown       }
271544be873SJed Brown     }
272544be873SJed Brown   } else {
2732b730f8bSJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
274544be873SJed Brown   }
275544be873SJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
276544be873SJed Brown   free(rows);
277544be873SJed Brown   free(cols);
278544be873SJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
279544be873SJed Brown   PetscFunctionReturn(0);
280544be873SJed Brown }
281544be873SJed Brown 
2822b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
283544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
284544be873SJed Brown   const PetscScalar *values;
285544be873SJed Brown   MatType            mat_type;
286544be873SJed Brown 
287544be873SJed Brown   PetscFunctionBeginUser;
288544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
2892b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
290544be873SJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
2912b730f8bSJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
292544be873SJed Brown   } else {
293544be873SJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
294544be873SJed Brown   }
295544be873SJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
296544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
297544be873SJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
298544be873SJed Brown   PetscFunctionReturn(0);
299544be873SJed Brown }
300544be873SJed Brown 
3012b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
302e334ad8fSJed Brown   User      user = *(User *)user_data;
303d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
304e334ad8fSJed Brown   PetscFunctionBeginUser;
3052b730f8bSJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorContextSetDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
306d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
307e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3082b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
309e334ad8fSJed Brown   if (!user->matrices_set_up) {
310e334ad8fSJed Brown     if (J_is_shell) {
311e334ad8fSJed Brown       PetscCall(MatShellSetContext(J, user));
3122b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3132b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
314e334ad8fSJed Brown       PetscCall(MatSetUp(J));
315e334ad8fSJed Brown     }
316e334ad8fSJed Brown     if (!J_pre_is_shell) {
3172b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
318544be873SJed Brown     }
319d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
320544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
321544be873SJed Brown     }
322e334ad8fSJed Brown     user->matrices_set_up = true;
323e334ad8fSJed Brown   }
324e334ad8fSJed Brown   if (!J_pre_is_shell) {
3252b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
326e334ad8fSJed Brown   }
327d6bf345cSJed Brown   if (user->coo_values_amat) {
328d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
329d6bf345cSJed Brown   } else if (J_is_mffd) {
330d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
331d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
332d6bf345cSJed Brown   }
333e334ad8fSJed Brown   PetscFunctionReturn(0);
334e334ad8fSJed Brown }
335e334ad8fSJed Brown 
3362b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
33777841947SLeila Ghaffari   Vec         Q_loc;
33877841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
33977841947SLeila Ghaffari   PetscViewer viewer;
34077841947SLeila Ghaffari   PetscFunctionBeginUser;
34177841947SLeila Ghaffari 
34237cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
34377841947SLeila Ghaffari     // Set up output
344f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
345f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
346f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
347f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
34877841947SLeila Ghaffari 
34977841947SLeila Ghaffari     // Output
35037cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
351f14660a4SJames Wright 
3522b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
353f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
354f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
35577841947SLeila Ghaffari     if (user->dm_viz) {
35677841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
35777841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
35877841947SLeila Ghaffari       PetscViewer viewer_refined;
35977841947SLeila Ghaffari 
360f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
361f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
362f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
363f14660a4SJames Wright 
364f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
365f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3662b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
367f14660a4SJames Wright 
36837cbb16aSJed Brown       PetscCall(
36937cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
370f14660a4SJames Wright 
3712b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
372f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
373f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
374f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
375f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
37677841947SLeila Ghaffari     }
377f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
37837cbb16aSJed Brown   }
37977841947SLeila Ghaffari 
38077841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
38169293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
38237cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
38369293791SJames Wright   } else {
3842b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
38569293791SJames Wright   }
3862b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
387f14660a4SJames Wright 
3884de8550aSJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
3894de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
3904de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3914de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
3924de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
393f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
394f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
395f14660a4SJames Wright   PetscFunctionReturn(0);
396f14660a4SJames Wright }
397f14660a4SJames Wright 
398f14660a4SJames Wright // User provided TS Monitor
3992b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
400f14660a4SJames Wright   User user = ctx;
401f14660a4SJames Wright   PetscFunctionBeginUser;
402f14660a4SJames Wright 
40337cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
404894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
405894de27cSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
406894de27cSJames Wright     PetscFunctionReturn(0);
407f14660a4SJames Wright 
408f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
40977841947SLeila Ghaffari 
41077841947SLeila Ghaffari   PetscFunctionReturn(0);
41177841947SLeila Ghaffari }
41277841947SLeila Ghaffari 
41377841947SLeila Ghaffari // TS: Create, setup, and solve
4142b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
41577841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
41677841947SLeila Ghaffari   TSAdapt     adapt;
41777841947SLeila Ghaffari   PetscScalar final_time;
41877841947SLeila Ghaffari   PetscFunctionBeginUser;
41977841947SLeila Ghaffari 
4202b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4212b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
42277841947SLeila Ghaffari   if (phys->implicit) {
4232b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
42477841947SLeila Ghaffari     if (user->op_ifunction) {
4252b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
42677841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4272b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
42877841947SLeila Ghaffari     }
429e334ad8fSJed Brown     if (user->op_ijacobian) {
4302b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
431544be873SJed Brown       if (app_ctx->amat_type) {
432544be873SJed Brown         Mat Pmat, Amat;
4332b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4342b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4352b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4362b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4372b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4382b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
439544be873SJed Brown       }
440e334ad8fSJed Brown     }
44177841947SLeila Ghaffari   } else {
4422b730f8bSJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4432b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4442b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4452b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
44677841947SLeila Ghaffari   }
4472b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4482b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
449cf7a0454SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4502b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4512b730f8bSJeremy L Thompson   if (app_ctx->test_mode) {
4522b730f8bSJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4532b730f8bSJeremy L Thompson   }
4542b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4552b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4562b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
4575e82a6e1SJeremy L Thompson   user->time = -1.0;  // require all BCs and ctx to be updated
4585e82a6e1SJeremy L Thompson   user->dt   = -1.0;
45977841947SLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
46077841947SLeila Ghaffari     if (!app_ctx->test_mode) {
4612b730f8bSJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
46277841947SLeila Ghaffari     }
46377841947SLeila Ghaffari   } else {  // continue from time of last output
46477841947SLeila Ghaffari     PetscInt    count;
46577841947SLeila Ghaffari     PetscViewer viewer;
4662b730f8bSJeremy L Thompson 
4674de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4682b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4694de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4702b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4714de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4724de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4734de8550aSJed Brown     }
4744de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
475d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
47677841947SLeila Ghaffari   }
47777841947SLeila Ghaffari   if (!app_ctx->test_mode) {
4782b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
47977841947SLeila Ghaffari   }
48077841947SLeila Ghaffari 
48177841947SLeila Ghaffari   // Solve
482d8e0aecdSJed Brown   PetscReal start_time;
483d8e0aecdSJed Brown   PetscInt  start_step;
4842b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
485d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
4867b39487dSJeremy L Thompson 
4877b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
4887b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
489d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
4907b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
4917b39487dSJeremy L Thompson     // LCOV_EXCL_START
4927b39487dSJeremy L Thompson     SNES      snes;
4937b39487dSJeremy L Thompson     Vec       Q_preload;
4947b39487dSJeremy L Thompson     PetscReal rtol;
4957b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
4967b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
4977b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
4987b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
4992b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5007b39487dSJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
5017b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
5022b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5037b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
5047b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
5057b39487dSJeremy L Thompson   } else {
5062b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5072b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
5087b39487dSJeremy L Thompson   }
5097b39487dSJeremy L Thompson   PetscPreLoadEnd();
5107b39487dSJeremy L Thompson 
5117b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
51277841947SLeila Ghaffari   *f_time = final_time;
5137b39487dSJeremy L Thompson 
51477841947SLeila Ghaffari   if (!app_ctx->test_mode) {
51537cbb16aSJed Brown     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
516f14660a4SJames Wright       PetscInt step_no;
517f14660a4SJames Wright       PetscCall(TSGetStepNumber(*ts, &step_no));
518f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
519f14660a4SJames Wright     }
520f14660a4SJames Wright 
5217b39487dSJeremy L Thompson     PetscLogEvent stage_id;
5227b39487dSJeremy L Thompson     PetscStageLog stage_log;
5237b39487dSJeremy L Thompson 
5247b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
5257b39487dSJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
5262b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
52777841947SLeila Ghaffari   }
52877841947SLeila Ghaffari   PetscFunctionReturn(0);
52977841947SLeila Ghaffari }
530