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 #include "../qfunctions/mass.h" 13a515125bSLeila Ghaffari 14a515125bSLeila Ghaffari // Compute mass matrix for explicit scheme 152b916ea7SJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) { 16a515125bSLeila Ghaffari Vec M_loc; 17a515125bSLeila Ghaffari CeedQFunction qf_mass; 18a515125bSLeila Ghaffari CeedOperator op_mass; 19a515125bSLeila Ghaffari CeedVector m_ceed, ones_vec; 20a515125bSLeila Ghaffari CeedInt num_comp_q, q_data_size; 21a515125bSLeila Ghaffari PetscFunctionBeginUser; 22a515125bSLeila Ghaffari 23a515125bSLeila Ghaffari // CEED Restriction 24a515125bSLeila Ghaffari CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q); 25a515125bSLeila Ghaffari CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size); 26a515125bSLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL); 27a515125bSLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL); 28a515125bSLeila Ghaffari CeedVectorSetValue(ones_vec, 1.0); 29a515125bSLeila Ghaffari 30a515125bSLeila Ghaffari // CEED QFunction 31a515125bSLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 32a515125bSLeila Ghaffari CeedQFunctionAddInput(qf_mass, "q", num_comp_q, CEED_EVAL_INTERP); 333b1e209bSJeremy L Thompson CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE); 34a515125bSLeila Ghaffari CeedQFunctionAddOutput(qf_mass, "v", num_comp_q, CEED_EVAL_INTERP); 35a515125bSLeila Ghaffari 36a515125bSLeila Ghaffari // CEED Operator 37a515125bSLeila Ghaffari CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 382b916ea7SJeremy L Thompson CeedOperatorSetField(op_mass, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 392b916ea7SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 402b916ea7SJeremy L Thompson CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 41a515125bSLeila Ghaffari 42a515125bSLeila Ghaffari // Place PETSc vector in CEED vector 43a515125bSLeila Ghaffari CeedScalar *m; 44a515125bSLeila Ghaffari PetscMemType m_mem_type; 452b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &M_loc)); 462b916ea7SJeremy L Thompson PetscCall(VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type)); 47a515125bSLeila Ghaffari CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m); 48a515125bSLeila Ghaffari 49a515125bSLeila Ghaffari // Apply CEED Operator 50a515125bSLeila Ghaffari CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE); 51a515125bSLeila Ghaffari 52a515125bSLeila Ghaffari // Restore vectors 53a515125bSLeila Ghaffari CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL); 542b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m)); 55a515125bSLeila Ghaffari 56a515125bSLeila Ghaffari // Local-to-Global 572b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(M)); 582b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M)); 592b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &M_loc)); 60a515125bSLeila Ghaffari 61a515125bSLeila Ghaffari // Invert diagonally lumped mass vector for RHS function 622b916ea7SJeremy L Thompson PetscCall(VecReciprocal(M)); 63a515125bSLeila Ghaffari 64a515125bSLeila Ghaffari // Cleanup 65a515125bSLeila Ghaffari CeedVectorDestroy(&ones_vec); 66a515125bSLeila Ghaffari CeedVectorDestroy(&m_ceed); 67a515125bSLeila Ghaffari CeedQFunctionDestroy(&qf_mass); 68a515125bSLeila Ghaffari CeedOperatorDestroy(&op_mass); 69a515125bSLeila Ghaffari 70a515125bSLeila Ghaffari PetscFunctionReturn(0); 71a515125bSLeila Ghaffari } 72a515125bSLeila Ghaffari 73a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 74a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 75a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 76a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 77a515125bSLeila Ghaffari User user = *(User *)user_data; 78a515125bSLeila Ghaffari PetscScalar *q, *g; 79e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, G_loc; 80a515125bSLeila Ghaffari PetscMemType q_mem_type, g_mem_type; 81a515125bSLeila Ghaffari PetscFunctionBeginUser; 82a515125bSLeila Ghaffari 83e2f84137SJeremy L Thompson // Get local vector 842b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(user->dm, &G_loc)); 85e2f84137SJeremy L Thompson 86e2f84137SJeremy L Thompson // Update time dependent data 87e2f84137SJeremy L Thompson if (user->time != t) { 882b916ea7SJeremy L Thompson PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 89e2f84137SJeremy L Thompson if (user->phys->solution_time_label) { 9092ada588SJames Wright CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t); 91e2f84137SJeremy L Thompson } 92e2f84137SJeremy L Thompson user->time = t; 93e2f84137SJeremy L Thompson } 94bb8a0c61SJames Wright if (user->phys->timestep_size_label) { 95bb8a0c61SJames Wright PetscScalar dt; 962b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 97e2f84137SJeremy L Thompson if (user->dt != dt) { 982b916ea7SJeremy L Thompson CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label, &dt); 99e2f84137SJeremy L Thompson user->dt = dt; 100e2f84137SJeremy L Thompson } 101bb8a0c61SJames Wright } 102a515125bSLeila Ghaffari 103a515125bSLeila Ghaffari // Global-to-local 1042b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 105a515125bSLeila Ghaffari 106a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 1072b916ea7SJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type)); 1082b916ea7SJeremy L Thompson PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type)); 109a515125bSLeila Ghaffari CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q); 110a515125bSLeila Ghaffari CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 111a515125bSLeila Ghaffari 112a515125bSLeila Ghaffari // Apply CEED operator 1132b916ea7SJeremy L Thompson CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 114a515125bSLeila Ghaffari 115a515125bSLeila Ghaffari // Restore vectors 116a515125bSLeila Ghaffari CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 117a515125bSLeila Ghaffari CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 1182b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q)); 1192b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(G_loc, &g)); 120a515125bSLeila Ghaffari 121a515125bSLeila Ghaffari // Local-to-Global 1222b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1232b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 124a515125bSLeila Ghaffari 125a515125bSLeila Ghaffari // Inverse of the lumped mass matrix (M is Minv) 1262b916ea7SJeremy L Thompson PetscCall(VecPointwiseMult(G, G, user->M)); 127a515125bSLeila Ghaffari 128a515125bSLeila Ghaffari // Restore vectors 1292b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(user->dm, &G_loc)); 130a515125bSLeila Ghaffari 131a515125bSLeila Ghaffari PetscFunctionReturn(0); 132a515125bSLeila Ghaffari } 133a515125bSLeila Ghaffari 134a515125bSLeila Ghaffari // Implicit time-stepper function setup 1352b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 136a515125bSLeila Ghaffari User user = *(User *)user_data; 137a515125bSLeila Ghaffari const PetscScalar *q, *q_dot; 138a515125bSLeila Ghaffari PetscScalar *g; 139e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 140a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 141a515125bSLeila Ghaffari PetscFunctionBeginUser; 142a515125bSLeila Ghaffari 143e2f84137SJeremy L Thompson // Get local vectors 1442b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(user->dm, &G_loc)); 145e2f84137SJeremy L Thompson 146e2f84137SJeremy L Thompson // Update time dependent data 147e2f84137SJeremy L Thompson if (user->time != t) { 1482b916ea7SJeremy L Thompson PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 149e2f84137SJeremy L Thompson if (user->phys->solution_time_label) { 1502b916ea7SJeremy L Thompson CeedOperatorContextSetDouble(user->op_ifunction, user->phys->solution_time_label, &t); 151e2f84137SJeremy L Thompson } 152e2f84137SJeremy L Thompson user->time = t; 153e2f84137SJeremy L Thompson } 154bb8a0c61SJames Wright if (user->phys->timestep_size_label) { 155bb8a0c61SJames Wright PetscScalar dt; 1562b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 157e2f84137SJeremy L Thompson if (user->dt != dt) { 1582b916ea7SJeremy L Thompson CeedOperatorContextSetDouble(user->op_ifunction, user->phys->timestep_size_label, &dt); 159e2f84137SJeremy L Thompson user->dt = dt; 160e2f84137SJeremy L Thompson } 161bb8a0c61SJames Wright } 162a515125bSLeila Ghaffari 163a515125bSLeila Ghaffari // Global-to-local 1642b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 1652b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 166a515125bSLeila Ghaffari 167a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 1682b916ea7SJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type)); 1692b916ea7SJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type)); 1702b916ea7SJeremy L Thompson PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type)); 1712b916ea7SJeremy L Thompson CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q); 1722b916ea7SJeremy L Thompson CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), CEED_USE_POINTER, (PetscScalar *)q_dot); 173a515125bSLeila Ghaffari CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 174a515125bSLeila Ghaffari 175a515125bSLeila Ghaffari // Apply CEED operator 1762b916ea7SJeremy L Thompson CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 177a515125bSLeila Ghaffari 178a515125bSLeila Ghaffari // Restore vectors 179a515125bSLeila Ghaffari CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 180a515125bSLeila Ghaffari CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL); 181a515125bSLeila Ghaffari CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 1822b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q)); 1832b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot)); 1842b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(G_loc, &g)); 185a515125bSLeila Ghaffari 186a515125bSLeila Ghaffari // Local-to-Global 1872b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1882b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 189a515125bSLeila Ghaffari 190a515125bSLeila Ghaffari // Restore vectors 1912b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(user->dm, &G_loc)); 192a515125bSLeila Ghaffari 193a515125bSLeila Ghaffari PetscFunctionReturn(0); 194a515125bSLeila Ghaffari } 195a515125bSLeila Ghaffari 196f0b65372SJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) { 197f0b65372SJed Brown User user; 198f0b65372SJed Brown const PetscScalar *q; 199f0b65372SJed Brown PetscScalar *g; 200f0b65372SJed Brown PetscMemType q_mem_type, g_mem_type; 201f0b65372SJed Brown PetscFunctionBeginUser; 2022b916ea7SJeremy L Thompson PetscCall(MatShellGetContext(J, &user)); 203e2f84137SJeremy L Thompson Vec Q_loc = user->Q_dot_loc, // Note - Q_dot_loc has zero BCs 204e2f84137SJeremy L Thompson G_loc; 205e2f84137SJeremy L Thompson 206f0b65372SJed Brown // Get local vectors 2072b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(user->dm, &G_loc)); 208f0b65372SJed Brown 209f0b65372SJed Brown // Global-to-local 2102b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 211f0b65372SJed Brown 212f0b65372SJed Brown // Place PETSc vectors in CEED vectors 2132b916ea7SJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type)); 2142b916ea7SJeremy L Thompson PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type)); 2152b916ea7SJeremy L Thompson CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q); 216f0b65372SJed Brown CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 217f0b65372SJed Brown 218f0b65372SJed Brown // Apply CEED operator 2192b916ea7SJeremy L Thompson CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 220f0b65372SJed Brown 221f0b65372SJed Brown // Restore vectors 222f0b65372SJed Brown CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 223f0b65372SJed Brown CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 2242b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q)); 2252b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(G_loc, &g)); 226f0b65372SJed Brown 227f0b65372SJed Brown // Local-to-Global 2282b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 2292b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 230f0b65372SJed Brown 231f0b65372SJed Brown // Restore vectors 2322b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(user->dm, &G_loc)); 233f0b65372SJed Brown PetscFunctionReturn(0); 234f0b65372SJed Brown } 235f0b65372SJed Brown 236f0b65372SJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) { 237f0b65372SJed Brown User user; 238f0b65372SJed Brown Vec D_loc; 239f0b65372SJed Brown PetscScalar *d; 240f0b65372SJed Brown PetscMemType mem_type; 241f0b65372SJed Brown 242f0b65372SJed Brown PetscFunctionBeginUser; 243f0b65372SJed Brown MatShellGetContext(A, &user); 244f0b65372SJed Brown PetscCall(DMGetLocalVector(user->dm, &D_loc)); 245f0b65372SJed Brown PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type)); 246f0b65372SJed Brown CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d); 2472b916ea7SJeremy L Thompson CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE); 248f0b65372SJed Brown CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL); 249f0b65372SJed Brown PetscCall(VecRestoreArrayAndMemType(D_loc, &d)); 250f0b65372SJed Brown PetscCall(VecZeroEntries(D)); 251f0b65372SJed Brown PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D)); 252f0b65372SJed Brown PetscCall(DMRestoreLocalVector(user->dm, &D_loc)); 253f0b65372SJed Brown VecViewFromOptions(D, NULL, "-diag_vec_view"); 254f0b65372SJed Brown PetscFunctionReturn(0); 255f0b65372SJed Brown } 256f0b65372SJed Brown 2572b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) { 258b107fddaSJed Brown PetscCount ncoo; 259b107fddaSJed Brown PetscInt *rows, *cols; 260b107fddaSJed Brown 261b107fddaSJed Brown PetscFunctionBeginUser; 262b107fddaSJed Brown if (pbdiagonal) { 263b107fddaSJed Brown CeedSize l_size; 264b107fddaSJed Brown CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL); 265b107fddaSJed Brown ncoo = l_size * 5; 266b107fddaSJed Brown rows = malloc(ncoo * sizeof(rows[0])); 267b107fddaSJed Brown cols = malloc(ncoo * sizeof(cols[0])); 268b107fddaSJed Brown for (PetscCount n = 0; n < l_size / 5; n++) { 269b107fddaSJed Brown for (PetscInt i = 0; i < 5; i++) { 270b107fddaSJed Brown for (PetscInt j = 0; j < 5; j++) { 271b107fddaSJed Brown rows[(n * 5 + i) * 5 + j] = n * 5 + i; 272b107fddaSJed Brown cols[(n * 5 + i) * 5 + j] = n * 5 + j; 273b107fddaSJed Brown } 274b107fddaSJed Brown } 275b107fddaSJed Brown } 276b107fddaSJed Brown } else { 2772b916ea7SJeremy L Thompson PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols)); 278b107fddaSJed Brown } 279b107fddaSJed Brown PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols)); 280b107fddaSJed Brown free(rows); 281b107fddaSJed Brown free(cols); 282b107fddaSJed Brown CeedVectorCreate(user->ceed, ncoo, coo_values); 283b107fddaSJed Brown PetscFunctionReturn(0); 284b107fddaSJed Brown } 285b107fddaSJed Brown 2862b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) { 287b107fddaSJed Brown CeedMemType mem_type = CEED_MEM_HOST; 288b107fddaSJed Brown const PetscScalar *values; 289b107fddaSJed Brown MatType mat_type; 290b107fddaSJed Brown 291b107fddaSJed Brown PetscFunctionBeginUser; 292b107fddaSJed Brown PetscCall(MatGetType(J, &mat_type)); 2932b916ea7SJeremy L Thompson if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 294b107fddaSJed Brown if (user->app_ctx->pmat_pbdiagonal) { 2952b916ea7SJeremy L Thompson CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE); 296b107fddaSJed Brown } else { 297b107fddaSJed Brown CeedOperatorLinearAssemble(user->op_ijacobian, coo_values); 298b107fddaSJed Brown } 299b107fddaSJed Brown CeedVectorGetArrayRead(coo_values, mem_type, &values); 300b107fddaSJed Brown PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES)); 301b107fddaSJed Brown CeedVectorRestoreArrayRead(coo_values, &values); 302b107fddaSJed Brown PetscFunctionReturn(0); 303b107fddaSJed Brown } 304b107fddaSJed Brown 3052b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 306f0b65372SJed Brown User user = *(User *)user_data; 30704855949SJed Brown PetscBool J_is_shell, J_is_mffd, J_pre_is_shell; 308f0b65372SJed Brown PetscFunctionBeginUser; 3092b916ea7SJeremy L Thompson if (user->phys->ijacobian_time_shift_label) CeedOperatorContextSetDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift); 31004855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 311f0b65372SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell)); 3122b916ea7SJeremy L Thompson PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell)); 313f0b65372SJed Brown if (!user->matrices_set_up) { 314f0b65372SJed Brown if (J_is_shell) { 315f0b65372SJed Brown PetscCall(MatShellSetContext(J, user)); 3162b916ea7SJeremy L Thompson PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian)); 3172b916ea7SJeremy L Thompson PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian)); 318f0b65372SJed Brown PetscCall(MatSetUp(J)); 319f0b65372SJed Brown } 320f0b65372SJed Brown if (!J_pre_is_shell) { 3212b916ea7SJeremy L Thompson PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat)); 322b107fddaSJed Brown } 32304855949SJed Brown if (J != J_pre && !J_is_shell && !J_is_mffd) { 324b107fddaSJed Brown PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat)); 325b107fddaSJed Brown } 326f0b65372SJed Brown user->matrices_set_up = true; 327f0b65372SJed Brown } 328f0b65372SJed Brown if (!J_pre_is_shell) { 3292b916ea7SJeremy L Thompson PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat)); 330f0b65372SJed Brown } 33104855949SJed Brown if (user->coo_values_amat) { 33204855949SJed Brown PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat)); 33304855949SJed Brown } else if (J_is_mffd) { 33404855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 33504855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 33604855949SJed Brown } 337f0b65372SJed Brown PetscFunctionReturn(0); 338f0b65372SJed Brown } 339f0b65372SJed Brown 3402b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 341a515125bSLeila Ghaffari Vec Q_loc; 342a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 343a515125bSLeila Ghaffari PetscViewer viewer; 344a515125bSLeila Ghaffari PetscFunctionBeginUser; 345a515125bSLeila Ghaffari 346852e5969SJed Brown if (user->app_ctx->checkpoint_vtk) { 347a515125bSLeila Ghaffari // Set up output 3487538d537SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 3497538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 3507538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 3517538d537SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 352a515125bSLeila Ghaffari 353a515125bSLeila Ghaffari // Output 354852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 3557538d537SJames Wright 3562b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 3577538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 3587538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 359a515125bSLeila Ghaffari if (user->dm_viz) { 360a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 361a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 362a515125bSLeila Ghaffari PetscViewer viewer_refined; 363a515125bSLeila Ghaffari 3647538d537SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 3657538d537SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 3667538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 3677538d537SJames Wright 3687538d537SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 3697538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 3702b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 3717538d537SJames Wright 372852e5969SJed Brown PetscCall( 373852e5969SJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 3747538d537SJames Wright 3752b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 3767538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 3777538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 3787538d537SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 3797538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 380a515125bSLeila Ghaffari } 3817538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 382852e5969SJed Brown } 383a515125bSLeila Ghaffari 384a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 38591a36801SJames Wright if (user->app_ctx->add_stepnum2bin) { 386852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 38791a36801SJames Wright } else { 3882b916ea7SJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 38991a36801SJames Wright } 3902b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 3917538d537SJames Wright 392*9293eaa1SJed Brown PetscInt token = FLUIDS_FILE_TOKEN; 393*9293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT)); 394*9293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 395*9293eaa1SJed Brown time /= user->units->second; // Dimensionalize time back 396*9293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 3977538d537SJames Wright PetscCall(VecView(Q, viewer)); 3987538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 3997538d537SJames Wright PetscFunctionReturn(0); 4007538d537SJames Wright } 4017538d537SJames Wright 4027538d537SJames Wright // User provided TS Monitor 4032b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 4047538d537SJames Wright User user = ctx; 4057538d537SJames Wright PetscFunctionBeginUser; 4067538d537SJames Wright 407852e5969SJed Brown // Print every 'checkpoint_interval' steps 408c539088bSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 409c539088bSJames Wright (user->app_ctx->cont_steps == step_no && step_no != 0)) 410c539088bSJames Wright PetscFunctionReturn(0); 4117538d537SJames Wright 4127538d537SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 413a515125bSLeila Ghaffari 414a515125bSLeila Ghaffari PetscFunctionReturn(0); 415a515125bSLeila Ghaffari } 416a515125bSLeila Ghaffari 417a515125bSLeila Ghaffari // TS: Create, setup, and solve 4182b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 419a515125bSLeila Ghaffari MPI_Comm comm = user->comm; 420a515125bSLeila Ghaffari TSAdapt adapt; 421a515125bSLeila Ghaffari PetscScalar final_time; 422a515125bSLeila Ghaffari PetscFunctionBeginUser; 423a515125bSLeila Ghaffari 4242b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 4252b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 426a515125bSLeila Ghaffari if (phys->implicit) { 4272b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 428a515125bSLeila Ghaffari if (user->op_ifunction) { 4292b916ea7SJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 430a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 4312b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 432a515125bSLeila Ghaffari } 433f0b65372SJed Brown if (user->op_ijacobian) { 4342b916ea7SJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 435b107fddaSJed Brown if (app_ctx->amat_type) { 436b107fddaSJed Brown Mat Pmat, Amat; 4372b916ea7SJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Pmat)); 4382b916ea7SJeremy L Thompson PetscCall(DMSetMatType(dm, app_ctx->amat_type)); 4392b916ea7SJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Amat)); 4402b916ea7SJeremy L Thompson PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL)); 4412b916ea7SJeremy L Thompson PetscCall(MatDestroy(&Amat)); 4422b916ea7SJeremy L Thompson PetscCall(MatDestroy(&Pmat)); 443b107fddaSJed Brown } 444f0b65372SJed Brown } 445a515125bSLeila Ghaffari } else { 4462b916ea7SJeremy L Thompson if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 4472b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 4482b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 4492b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 450a515125bSLeila Ghaffari } 4512b916ea7SJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 4522b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 4532b916ea7SJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 4542b916ea7SJeremy L Thompson if (app_ctx->test_mode) { 4552b916ea7SJeremy L Thompson PetscCall(TSSetMaxSteps(*ts, 10)); 4562b916ea7SJeremy L Thompson } 4572b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 4582b916ea7SJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 4592b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 460e2f84137SJeremy L Thompson user->time = -1.0; // require all BCs and ctx to be updated 461e2f84137SJeremy L Thompson user->dt = -1.0; 462a515125bSLeila Ghaffari if (!app_ctx->cont_steps) { // print initial condition 463a515125bSLeila Ghaffari if (!app_ctx->test_mode) { 4642b916ea7SJeremy L Thompson PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user)); 465a515125bSLeila Ghaffari } 466a515125bSLeila Ghaffari } else { // continue from time of last output 467a515125bSLeila Ghaffari PetscInt count; 468a515125bSLeila Ghaffari PetscViewer viewer; 4692b916ea7SJeremy L Thompson 470*9293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 4712b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 472*9293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 4732b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 474*9293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 475*9293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 476*9293eaa1SJed Brown } 477*9293eaa1SJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 47874a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 479a515125bSLeila Ghaffari } 480a515125bSLeila Ghaffari if (!app_ctx->test_mode) { 4812b916ea7SJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 482a515125bSLeila Ghaffari } 483a515125bSLeila Ghaffari 484a515125bSLeila Ghaffari // Solve 48574a6f4ddSJed Brown PetscReal start_time; 48674a6f4ddSJed Brown PetscInt start_step; 4872b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 48874a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 48991982731SJeremy L Thompson 49091982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 49191982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 49274a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 49391982731SJeremy L Thompson if (PetscPreLoadingOn) { 49491982731SJeremy L Thompson // LCOV_EXCL_START 49591982731SJeremy L Thompson SNES snes; 49691982731SJeremy L Thompson Vec Q_preload; 49791982731SJeremy L Thompson PetscReal rtol; 49891982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 49991982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 50091982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 50191982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 5022b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 50391982731SJeremy L Thompson PetscCall(TSSetSolution(*ts, *Q)); 50491982731SJeremy L Thompson PetscCall(TSStep(*ts)); 5052b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 50691982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 50791982731SJeremy L Thompson // LCOV_EXCL_STOP 50891982731SJeremy L Thompson } else { 5092b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 5102b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 51191982731SJeremy L Thompson } 51291982731SJeremy L Thompson PetscPreLoadEnd(); 51391982731SJeremy L Thompson 51491982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 515a515125bSLeila Ghaffari *f_time = final_time; 51691982731SJeremy L Thompson 517a515125bSLeila Ghaffari if (!app_ctx->test_mode) { 518852e5969SJed Brown if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 5197538d537SJames Wright PetscInt step_no; 5207538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 5217538d537SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 5227538d537SJames Wright } 5237538d537SJames Wright 52491982731SJeremy L Thompson PetscLogEvent stage_id; 52591982731SJeremy L Thompson PetscStageLog stage_log; 52691982731SJeremy L Thompson 52791982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 52891982731SJeremy L Thompson PetscCall(PetscLogGetStageLog(&stage_log)); 5292b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time)); 530a515125bSLeila Ghaffari } 531a515125bSLeila Ghaffari PetscFunctionReturn(0); 532a515125bSLeila Ghaffari } 533