1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 #pragma once 4 5 #include <ceed.h> 6 #include <bc_definition.h> 7 #include <honee-file.h> 8 #include <log_events.h> 9 #include <mat-ceed.h> 10 #include <petsc-ceed-utils.h> 11 #include <petscts.h> 12 #include <stdbool.h> 13 14 #include <petsc_ops.h> 15 #include "../qfunctions/newtonian_types.h" 16 17 #if PETSC_VERSION_LT(3, 22, 0) 18 #error "PETSc v3.22 or later is required" 19 #endif 20 21 // ----------------------------------------------------------------------------- 22 // Enums 23 // ----------------------------------------------------------------------------- 24 25 // Euler - test cases 26 typedef enum { 27 EULER_TEST_ISENTROPIC_VORTEX = 0, 28 EULER_TEST_1 = 1, 29 EULER_TEST_2 = 2, 30 EULER_TEST_3 = 3, 31 EULER_TEST_4 = 4, 32 EULER_TEST_5 = 5, 33 } EulerTestType; 34 static const char *const EulerTestTypes[] = {"ISENTROPIC_VORTEX", "1", "2", "3", "4", "5", "EulerTestType", "EULER_TEST_", NULL}; 35 36 // Advection - Wind types 37 static const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL}; 38 39 // Advection - Initial Condition Types 40 static const char *const AdvDifICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "WAVE", 41 "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_", NULL}; 42 43 // Advection - Wave Types 44 static const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL}; 45 46 // Advection - Bubble Continuity Types 47 static const char *const AdvDifBubbleContinuityTypes[] = { 48 "SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_", NULL}; 49 50 // Stabilization methods 51 static const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 52 53 // Stabilization tau constants 54 static const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL}; 55 56 // Test mode type 57 typedef enum { 58 TESTTYPE_NONE = 0, 59 TESTTYPE_SOLVER = 1, 60 TESTTYPE_TURB_SPANSTATS = 2, 61 TESTTYPE_DIFF_FILTER = 3, 62 } TestType; 63 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "TestType", "TESTTYPE_", NULL}; 64 65 // Subgrid-Stress mode type 66 typedef enum { 67 SGS_MODEL_NONE = 0, 68 SGS_MODEL_DATA_DRIVEN = 1, 69 } SGSModelType; 70 static const char *const SGSModelTypes[] = {"NONE", "DATA_DRIVEN", "SGSModelType", "SGS_MODEL_", NULL}; 71 72 // Subgrid-Stress mode type 73 typedef enum { 74 SGS_MODEL_DD_FUSED = 0, 75 SGS_MODEL_DD_SEQENTIAL_CEED = 1, 76 SGS_MODEL_DD_SEQENTIAL_TORCH = 2, 77 } SGSModelDDImplementation; 78 static const char *const SGSModelDDImplementations[] = {"FUSED", "SEQUENTIAL_CEED", "SEQUENTIAL_TORCH", "SGSModelDDImplementation", "SGS_MODEL_DD_", 79 NULL}; 80 81 // Mesh transformation type 82 typedef enum { 83 MESH_TRANSFORM_NONE = 0, 84 MESH_TRANSFORM_PLATEMESH = 1, 85 } MeshTransformType; 86 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL}; 87 88 static const char *const DifferentialFilterDampingFunctions[] = { 89 "NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL}; 90 91 static const char *const DivDiffFluxProjectionMethods[] = {"NONE", "DIRECT", "INDIRECT", "DivDiffFluxProjectionMethod", "DIV_DIFF_FLUX_PROJ_", NULL}; 92 93 // ----------------------------------------------------------------------------- 94 // Structs 95 // ----------------------------------------------------------------------------- 96 // Structs declarations 97 typedef struct AppCtx_private *AppCtx; 98 typedef struct Honee_private *Honee; 99 typedef struct Units_private *Units; 100 typedef struct SimpleBC_private *SimpleBC; 101 typedef struct Physics_private *Physics; 102 typedef struct ProblemData_private *ProblemData; 103 104 // Application context from user command line options 105 struct AppCtx_private { 106 // libCEED arguments 107 char ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend 108 PetscInt degree; 109 PetscInt q_extra; 110 // Solver arguments 111 MatType amat_type; 112 // Post-processing arguments 113 PetscInt checkpoint_interval; 114 PetscInt viz_refine; 115 PetscInt cont_steps; 116 PetscReal cont_time; 117 char cont_file[PETSC_MAX_PATH_LEN]; 118 char output_dir[PETSC_MAX_PATH_LEN]; 119 PetscBool add_stepnum2bin; 120 PetscBool checkpoint_vtk; 121 // Problem type arguments 122 PetscFunctionList problems; 123 char problem_name[PETSC_MAX_PATH_LEN]; 124 // Test mode arguments 125 TestType test_type; 126 PetscScalar test_tol; 127 char test_file_path[PETSC_MAX_PATH_LEN]; 128 // Turbulent spanwise statistics 129 PetscBool turb_spanstats_enable; 130 PetscInt turb_spanstats_collect_interval; 131 PetscInt turb_spanstats_viewer_interval; 132 PetscViewer turb_spanstats_viewer; 133 PetscViewerFormat turb_spanstats_viewer_format; 134 // Wall forces 135 struct { 136 PetscInt num_wall; 137 PetscInt *walls; 138 PetscViewer viewer; 139 PetscViewerFormat viewer_format; 140 PetscBool header_written; 141 } wall_forces; 142 // Subgrid Stress Model 143 SGSModelType sgs_model_type; 144 PetscBool sgs_train_enable; 145 // Differential Filtering 146 PetscBool diff_filter_monitor; 147 MeshTransformType mesh_transform_type; 148 // Divergence of Diffusive Flux Projection 149 DivDiffFluxProjectionMethod divFdiffproj_method; 150 }; 151 152 typedef struct { 153 DM dm; 154 PetscSF sf; // For communicating child data to parents 155 OperatorApplyContext op_stats_collect_ctx, op_proj_rhs_ctx; 156 PetscInt num_comp_stats; 157 Vec Child_Stats_loc, Parent_Stats_loc; 158 KSP ksp; // For the L^2 projection solve 159 CeedScalar span_width; // spanwise width of the child domain 160 PetscBool do_mms_test; 161 OperatorApplyContext mms_error_ctx; 162 CeedContextFieldLabel solution_time_label, previous_time_label; 163 } SpanStatsData; 164 165 typedef struct { 166 DM dm; 167 PetscInt num_comp; 168 OperatorApplyContext l2_rhs_ctx; 169 KSP ksp; 170 } *NodalProjectionData; 171 172 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData; 173 174 struct DivDiffFluxProjectionData_private { 175 PetscInt num_diff_flux_comps; 176 DivDiffFluxProjectionMethod method; 177 NodalProjectionData projection; 178 179 // CeedOperator Objects 180 CeedElemRestriction elem_restr_div_diff_flux; 181 CeedBasis basis_div_diff_flux; 182 CeedEvalMode eval_mode_div_diff_flux; 183 CeedVector div_diff_flux_ceed; 184 185 // Problem specific setup functions 186 PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *); 187 PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *); 188 189 // Only used for direct method: 190 Vec DivDiffFlux_loc; 191 PetscMemType DivDiffFlux_memtype; 192 PetscBool ceed_vec_has_array; 193 194 // Only used for indirect method: 195 OperatorApplyContext calc_div_diff_flux; 196 }; 197 198 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 199 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 200 typedef struct { 201 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 202 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 203 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 204 CeedVector sgs_nodal_ceed, grad_velo_ceed; 205 SgsDDNodalStressEval sgs_nodal_eval; 206 SgsDDNodalStressInference sgs_nodal_inference; 207 void *sgs_nodal_inference_ctx; 208 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 209 } *SgsDDData; 210 211 typedef struct { 212 DM dm_dd_training; 213 PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 214 PetscScalar filter_widths[16]; 215 OperatorApplyContext op_training_data_calc_ctx; 216 NodalProjectionData filtered_grad_velo_proj; 217 size_t training_data_array_dims[2]; 218 PetscBool overwrite_training_data; 219 } *SGS_DD_TrainingData; 220 221 typedef struct { 222 DM dm_filter; 223 PetscInt num_filtered_fields; 224 CeedInt *num_field_components; 225 PetscInt field_prim_state, field_velo_prod; 226 OperatorApplyContext op_rhs_ctx; 227 KSP ksp; 228 PetscObjectState X_loc_state; 229 PetscBool do_mms_test; 230 CeedContextFieldLabel filter_width_scaling_label; 231 } *DiffFilterData; 232 233 typedef struct { 234 void *client; 235 char rank_id_name[16]; 236 PetscInt collocated_database_num_ranks; 237 } *SmartSimData; 238 239 // PETSc user data 240 struct Honee_private { 241 MPI_Comm comm; 242 DM dm; 243 DM dm_viz; 244 Mat interp_viz; 245 Ceed ceed; 246 Units units; 247 Vec Q_loc, Q_dot_loc; 248 Physics phys; 249 AppCtx app_ctx; 250 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 251 CeedOperator op_ifunction; 252 Mat mat_ijacobian; 253 KSP mass_ksp; 254 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 255 CeedScalar time_bc_set; 256 SpanStatsData spanstats; 257 NodalProjectionData grad_velo_proj; 258 SgsDDData sgs_dd_data; 259 DiffFilterData diff_filter; 260 SmartSimData smartsim; 261 SGS_DD_TrainingData sgs_dd_train; 262 DivDiffFluxProjectionData diff_flux_proj; 263 264 // old CeedData 265 CeedVector x_coord; 266 CeedBasis basis_x, basis_q; 267 CeedElemRestriction elem_restr_x, elem_restr_q; 268 OperatorApplyContext op_ics_ctx; 269 }; 270 271 // Units 272 struct Units_private { 273 // fundamental units 274 PetscScalar meter; 275 PetscScalar kilogram; 276 PetscScalar second; 277 PetscScalar Kelvin; 278 // derived units 279 PetscScalar Pascal; 280 PetscScalar J_per_kg_K; 281 PetscScalar m_per_squared_s; 282 PetscScalar W_per_m_K; 283 PetscScalar Joule; 284 }; 285 286 // Boundary conditions 287 struct SimpleBC_private { 288 PetscInt num_inflow, num_outflow, num_freestream, num_slip; 289 PetscInt inflows[16], outflows[16], freestreams[16], slips[16]; 290 }; 291 292 // Struct that contains all enums and structs used for the physics of all problems 293 struct Physics_private { 294 PetscBool implicit; 295 StateVariable state_var; 296 CeedContextFieldLabel solution_time_label; 297 CeedContextFieldLabel stg_solution_time_label; 298 CeedContextFieldLabel timestep_size_label; 299 CeedContextFieldLabel ics_time_label; 300 }; 301 302 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx, SimpleBC bc); 303 304 typedef struct { 305 CeedQFunctionUser qf_func_ptr; // !< QFunction function pointer 306 const char *qf_loc; // !< Absolute path to QFunction source file 307 CeedQFunctionContext qfctx; // !< QFunctionContext to attach to QFunction 308 } ProblemQFunctionSpec; 309 310 // Problem specific data 311 struct ProblemData_private { 312 CeedInt jac_data_size_vol, jac_data_size_sur; 313 ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip, 314 apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian; 315 bool compute_exact_solution_error; 316 PetscBool set_bc_from_ics, use_strong_bc_ceed; 317 PetscCount num_bc_defs; 318 BCDefinition *bc_defs; 319 PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx); 320 PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *); 321 }; 322 323 extern int FreeContextPetsc(void *); 324 325 // ----------------------------------------------------------------------------- 326 // Set up problems 327 // ----------------------------------------------------------------------------- 328 // Set up function for each problem 329 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 330 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 331 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 332 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 333 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 334 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 335 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 336 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 337 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 338 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 339 340 // Print function for each problem 341 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx); 342 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx); 343 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx); 344 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx); 345 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx); 346 347 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts); 348 349 // ----------------------------------------------------------------------------- 350 // libCEED functions 351 // ----------------------------------------------------------------------------- 352 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 353 CeedElemRestriction *restriction); 354 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 355 CeedElemRestriction *restriction); 356 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 357 PetscInt q_data_size, CeedElemRestriction *restriction); 358 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 359 PetscInt q_data_size, CeedElemRestriction *restriction); 360 361 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis); 362 PetscErrorCode CreateCoordinateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, CeedBasis *basis); 363 PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, 364 CeedBasis *basis); 365 PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field, 366 CeedBasis *basis); 367 PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, char **face_label_name); 368 PetscErrorCode DMLabelCreateGlobalValueArray(DM dm, DMLabel label, PetscInt *num_values, PetscInt **value_array); 369 370 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc); 371 372 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 373 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 374 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size); 375 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 376 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 377 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size); 378 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size); 379 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord, 380 CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 381 PetscErrorCode QDataClearStoredData(); 382 // ----------------------------------------------------------------------------- 383 // Time-stepping functions 384 // ----------------------------------------------------------------------------- 385 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 386 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 387 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 388 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts); 389 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t); 390 391 // ----------------------------------------------------------------------------- 392 // Setup DM 393 // ----------------------------------------------------------------------------- 394 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 395 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys); 396 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys); 397 398 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 399 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 400 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 401 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 402 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 403 404 // ----------------------------------------------------------------------------- 405 // Process command line options 406 // ----------------------------------------------------------------------------- 407 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx); 408 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc); 409 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]); 410 411 // ----------------------------------------------------------------------------- 412 // Miscellaneous utility functions 413 // ----------------------------------------------------------------------------- 414 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 415 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 416 CeedVector *inv_multiplicity); 417 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time); 418 419 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 420 Vec grad_FVM); 421 422 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 423 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time); 424 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time); 425 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 426 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 427 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context); 428 429 // ----------------------------------------------------------------------------- 430 // Turbulence Statistics Collection Functions 431 // ----------------------------------------------------------------------------- 432 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, Honee honee, ProblemData problem); 433 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 434 PetscErrorCode TurbulenceStatisticsDestroy(Honee honee); 435 436 // ----------------------------------------------------------------------------- 437 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 438 // ----------------------------------------------------------------------------- 439 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem); 440 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 441 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc); 442 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input, 443 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 444 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 445 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 446 CeedVector *grid_aniso_vector); 447 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 448 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 449 450 // ----------------------------------------------------------------------------- 451 // Boundary Condition Related Functions 452 // ----------------------------------------------------------------------------- 453 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem, SimpleBC bc); 454 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 455 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 456 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 457 458 // ----------------------------------------------------------------------------- 459 // Differential Filtering Functions 460 // ----------------------------------------------------------------------------- 461 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, ProblemData problem); 462 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 463 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 464 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 465 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem); 466 467 // ----------------------------------------------------------------------------- 468 // SGS Data-Driven Training via SmartSim 469 // ----------------------------------------------------------------------------- 470 PetscErrorCode SmartSimSetup(Honee honee); 471 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 472 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem); 473 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 474 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 475 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 476 477 // ----------------------------------------------------------------------------- 478 // Divergence of Diffusive Flux Projection 479 // ----------------------------------------------------------------------------- 480 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj); 481 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 482 CeedVector *vector, CeedEvalMode *eval_mode); 483 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj); 484 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc); 485 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj); 486