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