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