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