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