1 // Copyright (c) 2017-2024, 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 #pragma once 8 9 #include <ceed-utils.h> 10 #include <ceed.h> 11 #include <mat-ceed.h> 12 #include <petscts.h> 13 #include <stdbool.h> 14 15 #include "./include/petsc_ops.h" 16 #include "qfunctions/newtonian_types.h" 17 18 #if PETSC_VERSION_LT(3, 21, 0) 19 #error "PETSc v3.21 or later is required" 20 #endif 21 22 // ----------------------------------------------------------------------------- 23 // Enums 24 // ----------------------------------------------------------------------------- 25 26 // Euler - test cases 27 typedef enum { 28 EULER_TEST_ISENTROPIC_VORTEX = 0, 29 EULER_TEST_1 = 1, 30 EULER_TEST_2 = 2, 31 EULER_TEST_3 = 3, 32 EULER_TEST_4 = 4, 33 EULER_TEST_5 = 5, 34 } EulerTestType; 35 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1", "test_2", "test_3", "test_4", "test_5", 36 "EulerTestType", "EULER_TEST_", NULL}; 37 38 // Advection - Wind types 39 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL}; 40 41 // Advection - Initial Condition Types 42 static const char *const AdvectionICTypes[] = {"sphere", "cylinder", "cosine_hill", "skew", "AdvectionICType", "ADVECTIONIC_", NULL}; 43 44 // Advection - Bubble Continuity Types 45 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "cosine", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL}; 46 47 // Stabilization methods 48 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 49 50 // Stabilization tau constants 51 static const char *const StabilizationTauTypes[] = {"Ctau", "AdvDiff_Shakib", "AdvDiff_Shakib_P", "StabilizationTauType", "STAB_TAU_", NULL}; 52 53 // Test mode type 54 typedef enum { 55 TESTTYPE_NONE = 0, 56 TESTTYPE_SOLVER = 1, 57 TESTTYPE_TURB_SPANSTATS = 2, 58 TESTTYPE_DIFF_FILTER = 3, 59 } TestType; 60 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "diff_filter", "TestType", "TESTTYPE_", NULL}; 61 62 // Subgrid-Stress mode type 63 typedef enum { 64 SGS_MODEL_NONE = 0, 65 SGS_MODEL_DATA_DRIVEN = 1, 66 } SGSModelType; 67 static const char *const SGSModelTypes[] = {"none", "data_driven", "SGSModelType", "SGS_MODEL_", NULL}; 68 69 // Mesh transformation type 70 typedef enum { 71 MESH_TRANSFORM_NONE = 0, 72 MESH_TRANSFORM_PLATEMESH = 1, 73 } MeshTransformType; 74 static const char *const MeshTransformTypes[] = {"none", "platemesh", "MeshTransformType", "MESH_TRANSFORM_", NULL}; 75 76 static const char *const DifferentialFilterDampingFunctions[] = { 77 "none", "van_driest", "mms", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL}; 78 79 // ----------------------------------------------------------------------------- 80 // Log Events 81 // ----------------------------------------------------------------------------- 82 extern PetscLogEvent FLUIDS_CeedOperatorApply; 83 extern PetscLogEvent FLUIDS_CeedOperatorAssemble; 84 extern PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 85 extern PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 86 extern PetscLogEvent FLUIDS_SmartRedis_Init; 87 extern PetscLogEvent FLUIDS_SmartRedis_Meta; 88 extern PetscLogEvent FLUIDS_SmartRedis_Train; 89 extern PetscLogEvent FLUIDS_TrainDataCompute; 90 extern PetscLogEvent FLUIDS_DifferentialFilter; 91 extern PetscLogEvent FLUIDS_VelocityGradientProjection; 92 PetscErrorCode RegisterLogEvents(); 93 94 // ----------------------------------------------------------------------------- 95 // Structs 96 // ----------------------------------------------------------------------------- 97 // Structs declarations 98 typedef struct AppCtx_private *AppCtx; 99 typedef struct CeedData_private *CeedData; 100 typedef struct User_private *User; 101 typedef struct Units_private *Units; 102 typedef struct SimpleBC_private *SimpleBC; 103 typedef struct Physics_private *Physics; 104 105 // Application context from user command line options 106 struct AppCtx_private { 107 // libCEED arguments 108 char ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend 109 PetscInt degree; 110 PetscInt q_extra; 111 // Solver arguments 112 MatType amat_type; 113 // Post-processing arguments 114 PetscInt checkpoint_interval; 115 PetscInt viz_refine; 116 PetscInt cont_steps; 117 PetscReal cont_time; 118 char cont_file[PETSC_MAX_PATH_LEN]; 119 char cont_time_file[PETSC_MAX_PATH_LEN]; 120 char output_dir[PETSC_MAX_PATH_LEN]; 121 PetscBool add_stepnum2bin; 122 PetscBool checkpoint_vtk; 123 // Problem type arguments 124 PetscFunctionList problems; 125 char problem_name[PETSC_MAX_PATH_LEN]; 126 // Test mode arguments 127 TestType test_type; 128 PetscScalar test_tol; 129 char test_file_path[PETSC_MAX_PATH_LEN]; 130 // Turbulent spanwise statistics 131 PetscBool turb_spanstats_enable; 132 PetscInt turb_spanstats_collect_interval; 133 PetscInt turb_spanstats_viewer_interval; 134 PetscViewer turb_spanstats_viewer; 135 PetscViewerFormat turb_spanstats_viewer_format; 136 // Wall forces 137 struct { 138 PetscInt num_wall; 139 PetscInt *walls; 140 PetscViewer viewer; 141 PetscViewerFormat viewer_format; 142 PetscBool header_written; 143 } wall_forces; 144 // Subgrid Stress Model 145 SGSModelType sgs_model_type; 146 PetscBool sgs_train_enable; 147 // Differential Filtering 148 PetscBool diff_filter_monitor; 149 MeshTransformType mesh_transform_type; 150 }; 151 152 // libCEED data struct 153 struct CeedData_private { 154 CeedVector x_coord, q_data; 155 CeedBasis basis_x, basis_q; 156 CeedElemRestriction elem_restr_x, elem_restr_q, elem_restr_qd_i; 157 OperatorApplyContext op_ics_ctx; 158 CeedQFunction qf_setup_sur; 159 }; 160 161 typedef struct { 162 DM dm; 163 PetscSF sf; // For communicating child data to parents 164 OperatorApplyContext op_stats_collect_ctx, op_proj_rhs_ctx; 165 PetscInt num_comp_stats; 166 Vec Child_Stats_loc, Parent_Stats_loc; 167 KSP ksp; // For the L^2 projection solve 168 CeedScalar span_width; // spanwise width of the child domain 169 PetscBool do_mms_test; 170 OperatorApplyContext mms_error_ctx; 171 CeedContextFieldLabel solution_time_label, previous_time_label; 172 } SpanStatsData; 173 174 typedef struct { 175 DM dm; 176 PetscInt num_comp; 177 OperatorApplyContext l2_rhs_ctx; 178 KSP ksp; 179 } *NodalProjectionData; 180 181 typedef PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 182 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 183 typedef struct { 184 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 185 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 186 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 187 CeedVector sgs_nodal_ceed, grad_velo_ceed; 188 SgsDDNodalStressEval sgs_nodal_eval; 189 SgsDDNodalStressInference sgs_nodal_inference; 190 void *sgs_nodal_inference_ctx; 191 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 192 } *SgsDDData; 193 194 typedef struct { 195 DM dm_dd_training; 196 PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 197 PetscScalar filter_widths[16]; 198 OperatorApplyContext op_training_data_calc_ctx; 199 NodalProjectionData filtered_grad_velo_proj; 200 size_t training_data_array_dims[2]; 201 PetscBool overwrite_training_data; 202 } *SGS_DD_TrainingData; 203 204 typedef struct { 205 DM dm_filter; 206 PetscInt num_filtered_fields; 207 CeedInt *num_field_components; 208 PetscInt field_prim_state, field_velo_prod; 209 OperatorApplyContext op_rhs_ctx; 210 KSP ksp; 211 PetscObjectState X_loc_state; 212 PetscBool do_mms_test; 213 CeedContextFieldLabel filter_width_scaling_label; 214 } *DiffFilterData; 215 216 typedef struct { 217 void *client; 218 char rank_id_name[16]; 219 PetscInt collocated_database_num_ranks; 220 } *SmartSimData; 221 222 // PETSc user data 223 struct User_private { 224 MPI_Comm comm; 225 DM dm; 226 DM dm_viz; 227 Mat interp_viz; 228 Ceed ceed; 229 Units units; 230 Vec Q_loc, Q_dot_loc; 231 Physics phys; 232 AppCtx app_ctx; 233 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 234 CeedOperator op_ifunction; 235 Mat mat_ijacobian; 236 KSP mass_ksp; 237 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 238 CeedScalar time_bc_set; 239 SpanStatsData spanstats; 240 NodalProjectionData grad_velo_proj; 241 SgsDDData sgs_dd_data; 242 DiffFilterData diff_filter; 243 SmartSimData smartsim; 244 SGS_DD_TrainingData sgs_dd_train; 245 }; 246 247 // Units 248 struct Units_private { 249 // fundamental units 250 PetscScalar meter; 251 PetscScalar kilogram; 252 PetscScalar second; 253 PetscScalar Kelvin; 254 // derived units 255 PetscScalar Pascal; 256 PetscScalar J_per_kg_K; 257 PetscScalar m_per_squared_s; 258 PetscScalar W_per_m_K; 259 PetscScalar Joule; 260 }; 261 262 // Boundary conditions 263 struct SimpleBC_private { 264 PetscInt num_wall, // Number of faces with wall BCs 265 wall_comps[5], // An array of constrained component numbers 266 num_comps, 267 num_symmetry[3], // Number of faces with symmetry BCs 268 num_inflow, num_outflow, num_freestream, num_slip; 269 PetscInt walls[16], symmetries[3][16], inflows[16], outflows[16], freestreams[16], slips[16]; 270 }; 271 272 // Struct that contains all enums and structs used for the physics of all problems 273 struct Physics_private { 274 PetscBool implicit; 275 StateVariable state_var; 276 CeedContextFieldLabel solution_time_label; 277 CeedContextFieldLabel stg_solution_time_label; 278 CeedContextFieldLabel timestep_size_label; 279 CeedContextFieldLabel ics_time_label; 280 CeedContextFieldLabel ijacobian_time_shift_label; 281 }; 282 283 typedef struct { 284 CeedQFunctionUser qfunction; 285 const char *qfunction_loc; 286 CeedQFunctionContext qfunction_context; 287 } ProblemQFunctionSpec; 288 289 // Problem specific data 290 typedef struct ProblemData_private *ProblemData; 291 struct ProblemData_private { 292 CeedInt dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur; 293 CeedScalar dm_scale; 294 ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, 295 apply_freestream, apply_slip, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian; 296 bool non_zero_time; 297 PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian; 298 PetscErrorCode (*print_info)(User, ProblemData, AppCtx); 299 PetscErrorCode (*create_mass_operator)(User, CeedOperator *); 300 }; 301 302 extern int FreeContextPetsc(void *); 303 304 // ----------------------------------------------------------------------------- 305 // Set up problems 306 // ----------------------------------------------------------------------------- 307 // Set up function for each problem 308 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 309 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 310 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 311 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 312 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 313 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 314 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 315 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 316 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 317 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 318 319 // Print function for each problem 320 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx); 321 322 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx); 323 324 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx); 325 326 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx); 327 328 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData problem, AppCtx app_ctx); 329 330 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, MPI_Comm comm); 331 332 // ----------------------------------------------------------------------------- 333 // libCEED functions 334 // ----------------------------------------------------------------------------- 335 // Utility function to create local CEED restriction 336 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field, 337 CeedElemRestriction *elem_restr); 338 339 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 340 CeedElemRestriction *restriction); 341 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 342 CeedElemRestriction *restriction); 343 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 344 PetscInt q_data_size, CeedElemRestriction *restriction); 345 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 346 PetscInt q_data_size, CeedElemRestriction *restriction); 347 348 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis); 349 350 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc); 351 352 // ----------------------------------------------------------------------------- 353 // Time-stepping functions 354 // ----------------------------------------------------------------------------- 355 // RHS (Explicit time-stepper) function setup 356 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 357 358 // Implicit time-stepper function setup 359 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 360 361 // User provided TS Monitor 362 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 363 364 // TS: Create, setup, and solve 365 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts); 366 367 // Update Boundary Values when time has changed 368 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t); 369 370 // ----------------------------------------------------------------------------- 371 // Setup DM 372 // ----------------------------------------------------------------------------- 373 // Create mesh 374 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 375 376 // Set up DM 377 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys); 378 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 379 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 380 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 381 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 382 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 383 384 // Refine DM for high-order viz 385 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys); 386 387 // ----------------------------------------------------------------------------- 388 // Process command line options 389 // ----------------------------------------------------------------------------- 390 // Register problems to be available on the command line 391 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx); 392 393 // Process general command line options 394 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc); 395 396 // ----------------------------------------------------------------------------- 397 // Miscellaneous utility functions 398 // ----------------------------------------------------------------------------- 399 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 400 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 401 CeedVector *inv_multiplicity); 402 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time); 403 404 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 405 Vec grad_FVM); 406 407 // Compare reference solution values with current test run for CI 408 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 409 410 // Get error for problems with exact solutions 411 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time); 412 413 // Post-processing 414 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time); 415 416 // -- Gather initial Q values in case of continuation of simulation 417 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q); 418 419 // Record boundary values from initial condition 420 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 421 422 // Versioning token for binary checkpoints 423 extern const PetscInt32 FLUIDS_FILE_TOKEN; // for backwards compatibility 424 extern const PetscInt32 FLUIDS_FILE_TOKEN_32; 425 extern const PetscInt32 FLUIDS_FILE_TOKEN_64; 426 427 // Create appropriate mass qfunction based on number of components N 428 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 429 430 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context); 431 432 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 433 FILE **fp); 434 435 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows); 436 437 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]); 438 439 // ----------------------------------------------------------------------------- 440 // Turbulence Statistics Collection Functions 441 // ----------------------------------------------------------------------------- 442 443 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 444 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 445 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data); 446 447 // ----------------------------------------------------------------------------- 448 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 449 // ----------------------------------------------------------------------------- 450 451 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 452 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 453 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc); 454 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input, 455 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 456 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 457 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 458 CeedVector *grid_aniso_vector); 459 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 460 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 461 462 // ----------------------------------------------------------------------------- 463 // Boundary Condition Related Functions 464 // ----------------------------------------------------------------------------- 465 466 // Setup StrongBCs that use QFunctions 467 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc); 468 469 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 470 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 471 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 472 473 // ----------------------------------------------------------------------------- 474 // Differential Filtering Functions 475 // ----------------------------------------------------------------------------- 476 477 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 478 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 479 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 480 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 481 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem); 482 483 // ----------------------------------------------------------------------------- 484 // SGS Data-Driven Training via SmartSim 485 // ----------------------------------------------------------------------------- 486 PetscErrorCode SmartSimSetup(User user); 487 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 488 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 489 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 490 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 491 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 492