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 }; 159 160 typedef struct { 161 DM dm; 162 PetscSF sf; // For communicating child data to parents 163 OperatorApplyContext op_stats_collect_ctx, op_proj_rhs_ctx; 164 PetscInt num_comp_stats; 165 Vec Child_Stats_loc, Parent_Stats_loc; 166 KSP ksp; // For the L^2 projection solve 167 CeedScalar span_width; // spanwise width of the child domain 168 PetscBool do_mms_test; 169 OperatorApplyContext mms_error_ctx; 170 CeedContextFieldLabel solution_time_label, previous_time_label; 171 } SpanStatsData; 172 173 typedef struct { 174 DM dm; 175 PetscInt num_comp; 176 OperatorApplyContext l2_rhs_ctx; 177 KSP ksp; 178 } *NodalProjectionData; 179 180 typedef PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 181 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 182 typedef struct { 183 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 184 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 185 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 186 CeedVector sgs_nodal_ceed, grad_velo_ceed; 187 SgsDDNodalStressEval sgs_nodal_eval; 188 SgsDDNodalStressInference sgs_nodal_inference; 189 void *sgs_nodal_inference_ctx; 190 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 191 } *SgsDDData; 192 193 typedef struct { 194 DM dm_dd_training; 195 PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 196 PetscScalar filter_widths[16]; 197 OperatorApplyContext op_training_data_calc_ctx; 198 NodalProjectionData filtered_grad_velo_proj; 199 size_t training_data_array_dims[2]; 200 PetscBool overwrite_training_data; 201 } *SGS_DD_TrainingData; 202 203 typedef struct { 204 DM dm_filter; 205 PetscInt num_filtered_fields; 206 CeedInt *num_field_components; 207 PetscInt field_prim_state, field_velo_prod; 208 OperatorApplyContext op_rhs_ctx; 209 KSP ksp; 210 PetscObjectState X_loc_state; 211 PetscBool do_mms_test; 212 CeedContextFieldLabel filter_width_scaling_label; 213 } *DiffFilterData; 214 215 typedef struct { 216 void *client; 217 char rank_id_name[16]; 218 PetscInt collocated_database_num_ranks; 219 } *SmartSimData; 220 221 // PETSc user data 222 struct User_private { 223 MPI_Comm comm; 224 DM dm; 225 DM dm_viz; 226 Mat interp_viz; 227 Ceed ceed; 228 Units units; 229 Vec Q_loc, Q_dot_loc; 230 Physics phys; 231 AppCtx app_ctx; 232 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 233 CeedOperator op_ifunction; 234 Mat mat_ijacobian; 235 KSP mass_ksp; 236 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 237 CeedScalar time_bc_set; 238 SpanStatsData spanstats; 239 NodalProjectionData grad_velo_proj; 240 SgsDDData sgs_dd_data; 241 DiffFilterData diff_filter; 242 SmartSimData smartsim; 243 SGS_DD_TrainingData sgs_dd_train; 244 }; 245 246 // Units 247 struct Units_private { 248 // fundamental units 249 PetscScalar meter; 250 PetscScalar kilogram; 251 PetscScalar second; 252 PetscScalar Kelvin; 253 // derived units 254 PetscScalar Pascal; 255 PetscScalar J_per_kg_K; 256 PetscScalar m_per_squared_s; 257 PetscScalar W_per_m_K; 258 PetscScalar Joule; 259 }; 260 261 // Boundary conditions 262 struct SimpleBC_private { 263 PetscInt num_wall, // Number of faces with wall BCs 264 wall_comps[5], // An array of constrained component numbers 265 num_comps, 266 num_symmetry[3], // Number of faces with symmetry BCs 267 num_inflow, num_outflow, num_freestream, num_slip; 268 PetscInt walls[16], symmetries[3][16], inflows[16], outflows[16], freestreams[16], slips[16]; 269 }; 270 271 // Struct that contains all enums and structs used for the physics of all problems 272 struct Physics_private { 273 PetscBool implicit; 274 StateVariable state_var; 275 CeedContextFieldLabel solution_time_label; 276 CeedContextFieldLabel stg_solution_time_label; 277 CeedContextFieldLabel timestep_size_label; 278 CeedContextFieldLabel ics_time_label; 279 CeedContextFieldLabel ijacobian_time_shift_label; 280 }; 281 282 typedef struct { 283 CeedQFunctionUser qfunction; 284 const char *qfunction_loc; 285 CeedQFunctionContext qfunction_context; 286 } ProblemQFunctionSpec; 287 288 // Problem specific data 289 typedef struct ProblemData_private *ProblemData; 290 struct ProblemData_private { 291 CeedInt dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur; 292 CeedScalar dm_scale; 293 ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip, 294 apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian; 295 bool non_zero_time; 296 PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian; 297 PetscErrorCode (*print_info)(User, ProblemData, AppCtx); 298 PetscErrorCode (*create_mass_operator)(User, CeedOperator *); 299 }; 300 301 extern int FreeContextPetsc(void *); 302 303 // ----------------------------------------------------------------------------- 304 // Set up problems 305 // ----------------------------------------------------------------------------- 306 // Set up function for each problem 307 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 308 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 309 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 310 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 311 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 312 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 313 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 314 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 315 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 316 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 317 318 // Print function for each problem 319 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx); 320 321 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx); 322 323 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx); 324 325 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx); 326 327 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData problem, AppCtx app_ctx); 328 329 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, MPI_Comm comm); 330 331 // ----------------------------------------------------------------------------- 332 // libCEED functions 333 // ----------------------------------------------------------------------------- 334 // Utility function to create local CEED restriction 335 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field, 336 CeedElemRestriction *elem_restr); 337 338 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 339 CeedElemRestriction *restriction); 340 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 341 CeedElemRestriction *restriction); 342 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 343 PetscInt q_data_size, CeedElemRestriction *restriction); 344 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 345 PetscInt q_data_size, CeedElemRestriction *restriction); 346 347 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis); 348 349 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc); 350 351 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 352 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 353 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size); 354 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 355 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 356 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size); 357 // ----------------------------------------------------------------------------- 358 // Time-stepping functions 359 // ----------------------------------------------------------------------------- 360 // RHS (Explicit time-stepper) function setup 361 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 362 363 // Implicit time-stepper function setup 364 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 365 366 // User provided TS Monitor 367 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 368 369 // TS: Create, setup, and solve 370 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts); 371 372 // Update Boundary Values when time has changed 373 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t); 374 375 // ----------------------------------------------------------------------------- 376 // Setup DM 377 // ----------------------------------------------------------------------------- 378 // Create mesh 379 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 380 381 // Set up DM 382 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys); 383 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 384 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 385 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 386 PetscErrorCode DMSetupByOrder_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 389 // Refine DM for high-order viz 390 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys); 391 392 // ----------------------------------------------------------------------------- 393 // Process command line options 394 // ----------------------------------------------------------------------------- 395 // Register problems to be available on the command line 396 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx); 397 398 // Process general command line options 399 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc); 400 401 // ----------------------------------------------------------------------------- 402 // Miscellaneous utility functions 403 // ----------------------------------------------------------------------------- 404 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 405 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 406 CeedVector *inv_multiplicity); 407 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time); 408 409 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 410 Vec grad_FVM); 411 412 // Compare reference solution values with current test run for CI 413 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 414 415 // Get error for problems with exact solutions 416 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time); 417 418 // Post-processing 419 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time); 420 421 // -- Gather initial Q values in case of continuation of simulation 422 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q); 423 424 // Record boundary values from initial condition 425 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 426 427 // Versioning token for binary checkpoints 428 extern const PetscInt32 FLUIDS_FILE_TOKEN; // for backwards compatibility 429 extern const PetscInt32 FLUIDS_FILE_TOKEN_32; 430 extern const PetscInt32 FLUIDS_FILE_TOKEN_64; 431 432 // Create appropriate mass qfunction based on number of components N 433 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 434 435 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context); 436 437 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 438 FILE **fp); 439 440 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows); 441 442 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]); 443 444 // ----------------------------------------------------------------------------- 445 // Turbulence Statistics Collection Functions 446 // ----------------------------------------------------------------------------- 447 448 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 449 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 450 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data); 451 452 // ----------------------------------------------------------------------------- 453 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 454 // ----------------------------------------------------------------------------- 455 456 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 457 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 458 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc); 459 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input, 460 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 461 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 462 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 463 CeedVector *grid_aniso_vector); 464 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 465 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 466 467 // ----------------------------------------------------------------------------- 468 // Boundary Condition Related Functions 469 // ----------------------------------------------------------------------------- 470 471 // Setup StrongBCs that use QFunctions 472 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc); 473 474 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 475 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 476 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 477 478 // ----------------------------------------------------------------------------- 479 // Differential Filtering Functions 480 // ----------------------------------------------------------------------------- 481 482 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 483 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 484 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 485 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 486 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem); 487 488 // ----------------------------------------------------------------------------- 489 // SGS Data-Driven Training via SmartSim 490 // ----------------------------------------------------------------------------- 491 PetscErrorCode SmartSimSetup(User user); 492 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 493 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 494 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 495 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 496 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 497