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 PetscReal initial_solution_time; 170 PetscInt initial_solution_step; 171 } *SpanStatsData; 172 173 typedef struct { 174 DM dm; 175 PetscInt num_comp; 176 OperatorApplyContext l2_rhs_ctx; 177 KSP ksp; 178 } *NodalProjectionData; 179 180 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData; 181 182 struct DivDiffFluxProjectionData_private { 183 PetscInt num_diff_flux_comps; 184 DivDiffFluxProjectionMethod method; 185 NodalProjectionData projection; 186 187 // CeedOperator Objects 188 CeedElemRestriction elem_restr_div_diff_flux; 189 CeedBasis basis_div_diff_flux; 190 CeedEvalMode eval_mode_div_diff_flux; 191 CeedVector div_diff_flux_ceed; 192 193 // Problem specific setup functions 194 PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *); 195 PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *); 196 197 // Only used for direct method: 198 Vec DivDiffFlux_loc; 199 PetscMemType DivDiffFlux_memtype; 200 PetscBool ceed_vec_has_array; 201 202 // Only used for indirect method: 203 OperatorApplyContext calc_div_diff_flux; 204 }; 205 206 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 207 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 208 typedef struct { 209 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 210 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 211 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 212 CeedVector sgs_nodal_ceed, grad_velo_ceed; 213 SgsDDNodalStressEval sgs_nodal_eval; 214 SgsDDNodalStressInference sgs_nodal_inference; 215 void *sgs_nodal_inference_ctx; 216 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 217 } *SgsDDData; 218 219 typedef struct { 220 DM dm_dd_training; 221 PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 222 PetscScalar filter_widths[16]; 223 OperatorApplyContext op_training_data_calc_ctx; 224 NodalProjectionData filtered_grad_velo_proj; 225 size_t training_data_array_dims[2]; 226 PetscBool overwrite_training_data; 227 } *SGS_DD_TrainingData; 228 229 typedef struct { 230 DM dm_filter; 231 PetscInt num_filtered_fields; 232 CeedInt *num_field_components; 233 PetscInt field_prim_state, field_velo_prod; 234 OperatorApplyContext op_rhs_ctx; 235 KSP ksp; 236 PetscObjectState X_loc_state; 237 PetscBool do_mms_test; 238 CeedContextFieldLabel filter_width_scaling_label; 239 } *DiffFilterData; 240 241 typedef struct { 242 void *client; 243 char rank_id_name[16]; 244 PetscInt collocated_database_num_ranks; 245 } *SmartSimData; 246 247 // PETSc user data 248 struct Honee_private { 249 MPI_Comm comm; 250 DM dm; 251 DM dm_viz; 252 Mat interp_viz; 253 Ceed ceed; 254 Units units; 255 Vec Q_loc, Q_dot_loc; 256 Physics phys; 257 AppCtx app_ctx; 258 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 259 CeedOperator op_ifunction; 260 Mat mat_ijacobian; 261 KSP mass_ksp; 262 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 263 CeedScalar time_bc_set; 264 SpanStatsData spanstats; 265 NodalProjectionData grad_velo_proj; 266 SgsDDData sgs_dd_data; 267 DiffFilterData diff_filter; 268 SmartSimData smartsim; 269 SGS_DD_TrainingData sgs_dd_train; 270 DivDiffFluxProjectionData diff_flux_proj; 271 272 // old CeedData 273 CeedVector x_coord; 274 CeedBasis basis_x, basis_q; 275 CeedElemRestriction elem_restr_x, elem_restr_q; 276 OperatorApplyContext op_ics_ctx; 277 278 PetscBool set_poststep; 279 time_t start_time; 280 time_t max_wall_time; 281 PetscInt max_wall_time_interval; 282 }; 283 284 // Units 285 struct Units_private { 286 // fundamental units 287 PetscScalar meter; 288 PetscScalar kilogram; 289 PetscScalar second; 290 PetscScalar Kelvin; 291 // derived units 292 PetscScalar Pascal; 293 PetscScalar J_per_kg_K; 294 PetscScalar m_per_squared_s; 295 PetscScalar W_per_m_K; 296 PetscScalar Joule; 297 }; 298 299 // Boundary conditions 300 struct SimpleBC_private { 301 PetscInt num_inflow, num_outflow, num_freestream, num_slip; 302 PetscInt inflows[16], outflows[16], freestreams[16], slips[16]; 303 }; 304 305 // Struct that contains all enums and structs used for the physics of all problems 306 struct Physics_private { 307 PetscBool implicit; 308 StateVariable state_var; 309 CeedContextFieldLabel solution_time_label; 310 CeedContextFieldLabel stg_solution_time_label; 311 CeedContextFieldLabel timestep_size_label; 312 CeedContextFieldLabel ics_time_label; 313 }; 314 315 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx, SimpleBC bc); 316 317 typedef struct { 318 CeedQFunctionUser qf_func_ptr; // !< QFunction function pointer 319 const char *qf_loc; // !< Absolute path to QFunction source file 320 CeedQFunctionContext qfctx; // !< QFunctionContext to attach to QFunction 321 } ProblemQFunctionSpec; 322 323 // Problem specific data 324 struct ProblemData_private { 325 CeedInt jac_data_size_vol, jac_data_size_sur; 326 ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip, 327 apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian; 328 bool compute_exact_solution_error; 329 PetscBool set_bc_from_ics, use_strong_bc_ceed; 330 PetscCount num_bc_defs; 331 BCDefinition *bc_defs; 332 PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx); 333 PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *); 334 }; 335 336 extern int FreeContextPetsc(void *); 337 338 // ----------------------------------------------------------------------------- 339 // Set up problems 340 // ----------------------------------------------------------------------------- 341 // Set up function for each problem 342 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 343 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 344 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 345 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 346 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 347 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 348 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 349 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 350 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 351 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 352 353 // Print function for each problem 354 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx); 355 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx); 356 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx); 357 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx); 358 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx); 359 360 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts); 361 362 // ----------------------------------------------------------------------------- 363 // libCEED functions 364 // ----------------------------------------------------------------------------- 365 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 366 CeedElemRestriction *restriction); 367 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 368 CeedElemRestriction *restriction); 369 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 370 PetscInt q_data_size, CeedElemRestriction *restriction); 371 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 372 PetscInt q_data_size, CeedElemRestriction *restriction); 373 374 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis); 375 PetscErrorCode CreateCoordinateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, CeedBasis *basis); 376 PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, 377 CeedBasis *basis); 378 PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field, 379 CeedBasis *basis); 380 PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, char **face_label_name); 381 PetscErrorCode DMLabelCreateGlobalValueArray(DM dm, DMLabel label, PetscInt *num_values, PetscInt **value_array); 382 383 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc); 384 385 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 386 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 387 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size); 388 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 389 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 390 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size); 391 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size); 392 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord, 393 CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 394 PetscErrorCode QDataClearStoredData(); 395 // ----------------------------------------------------------------------------- 396 // Time-stepping functions 397 // ----------------------------------------------------------------------------- 398 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 399 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 400 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 401 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts); 402 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t); 403 404 // ----------------------------------------------------------------------------- 405 // Setup DM 406 // ----------------------------------------------------------------------------- 407 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 408 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys); 409 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys); 410 411 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 412 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 413 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 414 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 415 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 416 417 // ----------------------------------------------------------------------------- 418 // Process command line options 419 // ----------------------------------------------------------------------------- 420 PetscErrorCode ProcessCommandLineOptions(Honee honee, SimpleBC bc); 421 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]); 422 423 // ----------------------------------------------------------------------------- 424 // Miscellaneous utility functions 425 // ----------------------------------------------------------------------------- 426 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 427 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 428 CeedVector *inv_multiplicity); 429 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time); 430 431 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 432 Vec grad_FVM); 433 434 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 435 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time); 436 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time); 437 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 438 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 439 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume); 440 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context); 441 442 // ----------------------------------------------------------------------------- 443 // Turbulence Statistics Collection Functions 444 // ----------------------------------------------------------------------------- 445 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, Honee honee, ProblemData problem); 446 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 447 PetscErrorCode TurbulenceStatisticsDestroy(SpanStatsData spanstats); 448 449 // ----------------------------------------------------------------------------- 450 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 451 // ----------------------------------------------------------------------------- 452 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem); 453 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 454 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc); 455 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input, 456 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 457 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 458 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 459 CeedVector *grid_aniso_vector); 460 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 461 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 462 463 // ----------------------------------------------------------------------------- 464 // Boundary Condition Related Functions 465 // ----------------------------------------------------------------------------- 466 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem, SimpleBC bc); 467 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 468 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 469 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 470 471 // ----------------------------------------------------------------------------- 472 // Differential Filtering Functions 473 // ----------------------------------------------------------------------------- 474 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, ProblemData problem); 475 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 476 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 477 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 478 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem); 479 480 // ----------------------------------------------------------------------------- 481 // SGS Data-Driven Training via SmartSim 482 // ----------------------------------------------------------------------------- 483 PetscErrorCode SmartSimSetup(Honee honee); 484 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 485 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem); 486 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 487 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 488 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 489 490 // ----------------------------------------------------------------------------- 491 // Divergence of Diffusive Flux Projection 492 // ----------------------------------------------------------------------------- 493 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj); 494 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 495 CeedVector *vector, CeedEvalMode *eval_mode); 496 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj); 497 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc); 498 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj); 499 500 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx); 501 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 502 503 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx); 504 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 505 506 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx); 507