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