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