Lines Matching +full:- +full:t

1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
4 // SPDX-License-Identifier: BSD-2-Clause
9 /// Time-stepping functions for Navier-Stokes example using PETSc
19 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { in UpdateBoundaryValues() argument
21 if (user->time_bc_set != t) { in UpdateBoundaryValues()
22 PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); in UpdateBoundaryValues()
23 user->time_bc_set = t; in UpdateBoundaryValues()
28 // RHS (Explicit time-stepper) function setup
29 // This is the RHS of the ODE, given as u_t = G(t,u)
31 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { in RHS_NS() argument
33 Ceed ceed = user->ceed; in RHS_NS()
35 Vec Q_loc = user->Q_loc; in RHS_NS()
40 PetscCall(UpdateBoundaryValues(user, Q_loc, t)); in RHS_NS()
41 …if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs… in RHS_NS()
43 …if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs… in RHS_NS()
45 PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); in RHS_NS()
47 PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); in RHS_NS()
50 PetscCall(KSPSolve(user->mass_ksp, G, G)); in RHS_NS()
52 PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); in RHS_NS()
88 reaction_force[w * dim + j] -= r[node].momentum[j]; in Surface_Forces_NS()
102 // Implicit time-stepper function setup
103 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { in IFunction_NS() argument
105 Ceed ceed = user->ceed; in IFunction_NS()
107 Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; in IFunction_NS()
112 PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); in IFunction_NS()
115 PetscCall(UpdateBoundaryValues(user, Q_loc, t)); in IFunction_NS()
116 …if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifu… in IFunction_NS()
118 …if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifu… in IFunction_NS()
120 // Global-to-local in IFunction_NS()
121 PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); in IFunction_NS()
122 PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); in IFunction_NS()
123 PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); in IFunction_NS()
124 PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); in IFunction_NS()
127 PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); in IFunction_NS()
128 PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); in IFunction_NS()
129 PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed)); in IFunction_NS()
134 …PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_R… in IFunction_NS()
139 PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); in IFunction_NS()
140 PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); in IFunction_NS()
141 PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc)); in IFunction_NS()
143 // Local-to-Global in IFunction_NS()
145 PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); in IFunction_NS()
148 PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); in IFunction_NS()
152 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J… in FormIJacobian_NS() argument
162 PetscCall(MatCeedSetContextReal(user->mat_ijacobian, "ijacobian time shift", shift)); in FormIJacobian_NS()
167 } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J)); in FormIJacobian_NS()
173 PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre)); in FormIJacobian_NS()
184 if (user->app_ctx->checkpoint_vtk) { in WriteOutput()
186 PetscCall(DMGetLocalVector(user->dm, &Q_loc)); in WriteOutput()
189 PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); in WriteOutput()
192 …scCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->… in WriteOutput()
197 if (user->dm_viz) { in WriteOutput()
202 PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); in WriteOutput()
203 PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); in WriteOutput()
206 PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); in WriteOutput()
208 PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); in WriteOutput()
210 …le_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->… in WriteOutput()
215 PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); in WriteOutput()
216 PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); in WriteOutput()
219 PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); in WriteOutput()
223 if (user->app_ctx->add_stepnum2bin) { in WriteOutput()
224 …ll(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ct… in WriteOutput()
226 …PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_d… in WriteOutput()
228 PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); in WriteOutput()
233 time /= user->units->second; // Dimensionalize time back in WriteOutput()
244 PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; in TSMonitor_WallForce()
245 const PetscInt *walls = user->app_ctx->wall_forces.walls; in TSMonitor_WallForce()
246 PetscViewer viewer = user->app_ctx->wall_forces.viewer; in TSMonitor_WallForce()
247 PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; in TSMonitor_WallForce()
253 PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); in TSMonitor_WallForce()
255 PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); in TSMonitor_WallForce()
256 PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); in TSMonitor_WallForce()
261 if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { in TSMonitor_WallForce()
263 user->app_ctx->wall_forces.header_written = PETSC_TRUE; in TSMonitor_WallForce()
287 …if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || in TSMonitor_NS()
288 (user->app_ctx->cont_steps == step_no && step_no != 0)) { in TSMonitor_NS()
298 MPI_Comm comm = user->comm; in TSSolve_NS()
306 if (phys->implicit) { in TSSolve_NS()
308 if (user->op_ifunction) { in TSSolve_NS()
313 if (user->mat_ijacobian) { in TSSolve_NS()
317 PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); in TSSolve_NS()
322 PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); in TSSolve_NS()
324 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); in TSSolve_NS()
325 PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); in TSSolve_NS()
327 PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); in TSSolve_NS()
329 if (user->mat_ijacobian) { in TSSolve_NS()
330 if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { in TSSolve_NS()
337 … PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); in TSSolve_NS()
338 PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL)); in TSSolve_NS()
343 user->time_bc_set = -1.0; // require all BCs be updated in TSSolve_NS()
344 if (app_ctx->cont_steps) { // continue from previous timestep data in TSSolve_NS()
348 if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time in TSSolve_NS()
349 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); in TSSolve_NS()
350 PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); in TSSolve_NS()
352 PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, in TSSolve_NS()
353 …"-continue step number not specified, but checkpoint file does not contain a step number (likely w… in TSSolve_NS()
355 PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); in TSSolve_NS()
356 PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); in TSSolve_NS()
358 if (app_ctx->test_type == TESTTYPE_NONE) { in TSSolve_NS()
361 if (app_ctx->wall_forces.viewer) { in TSSolve_NS()
364 if (app_ctx->turb_spanstats_enable) { in TSSolve_NS()
366 CeedScalar previous_time = app_ctx->cont_time * user->units->second; in TSSolve_NS()
367 PetscCallCeed(user->ceed, in TSSolve_NS()
368 …CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_ti… in TSSolve_NS()
370 …if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, … in TSSolve_NS()
372 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(user, user->phys, problem, *ts)); in TSSolve_NS()
379 PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view in TSSolve_NS()
407 if (app_ctx->test_type == TESTTYPE_NONE) { in TSSolve_NS()
410 if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { in TSSolve_NS()