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