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