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