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