1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 #pragma once 8 9 #include <ceed.h> 10 #include <log_events.h> 11 #include <mat-ceed.h> 12 #include <petsc-ceed-utils.h> 13 #include <petscts.h> 14 #include <stdbool.h> 15 16 #include "./include/petsc_ops.h" 17 #include "qfunctions/newtonian_types.h" 18 19 #if PETSC_VERSION_LT(3, 21, 0) 20 #error "PETSc v3.21 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 WindTypes[] = {"ROTATION", "TRANSLATION", "WindType", "WIND_", NULL}; 40 41 // Advection - Initial Condition Types 42 static const char *const AdvectionICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "AdvectionICType", "ADVECTIONIC_", NULL}; 43 44 // Advection - Bubble Continuity Types 45 static const char *const BubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL}; 46 47 // Stabilization methods 48 static const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 49 50 // Stabilization tau constants 51 static const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "ADVDIFF_SHAKIB_P", "StabilizationTauType", "STAB_TAU_", NULL}; 52 53 // Test mode type 54 typedef enum { 55 TESTTYPE_NONE = 0, 56 TESTTYPE_SOLVER = 1, 57 TESTTYPE_TURB_SPANSTATS = 2, 58 TESTTYPE_DIFF_FILTER = 3, 59 } TestType; 60 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "TestType", "TESTTYPE_", NULL}; 61 62 // Subgrid-Stress mode type 63 typedef enum { 64 SGS_MODEL_NONE = 0, 65 SGS_MODEL_DATA_DRIVEN = 1, 66 } SGSModelType; 67 static const char *const SGSModelTypes[] = {"NONE", "DATA_DRIVEN", "SGSModelType", "SGS_MODEL_", NULL}; 68 69 // Subgrid-Stress mode type 70 typedef enum { 71 SGS_MODEL_DD_FUSED = 0, 72 SGS_MODEL_DD_SEQENTIAL_CEED = 1, 73 SGS_MODEL_DD_SEQENTIAL_TORCH = 2, 74 } SGSModelDDImplementation; 75 static const char *const SGSModelDDImplementations[] = {"FUSED", "SEQUENTIAL_CEED", "SEQUENTIAL_TORCH", "SGSModelDDImplementation", "SGS_MODEL_DD_", 76 NULL}; 77 78 // Mesh transformation type 79 typedef enum { 80 MESH_TRANSFORM_NONE = 0, 81 MESH_TRANSFORM_PLATEMESH = 1, 82 } MeshTransformType; 83 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL}; 84 85 static const char *const DifferentialFilterDampingFunctions[] = { 86 "NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", 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 99 // Application context from user command line options 100 struct AppCtx_private { 101 // libCEED arguments 102 char ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend 103 PetscInt degree; 104 PetscInt q_extra; 105 // Solver arguments 106 MatType amat_type; 107 // Post-processing arguments 108 PetscInt checkpoint_interval; 109 PetscInt viz_refine; 110 PetscInt cont_steps; 111 PetscReal cont_time; 112 char cont_file[PETSC_MAX_PATH_LEN]; 113 char cont_time_file[PETSC_MAX_PATH_LEN]; 114 char output_dir[PETSC_MAX_PATH_LEN]; 115 PetscBool add_stepnum2bin; 116 PetscBool checkpoint_vtk; 117 // Problem type arguments 118 PetscFunctionList problems; 119 char problem_name[PETSC_MAX_PATH_LEN]; 120 // Test mode arguments 121 TestType test_type; 122 PetscScalar test_tol; 123 char test_file_path[PETSC_MAX_PATH_LEN]; 124 // Turbulent spanwise statistics 125 PetscBool turb_spanstats_enable; 126 PetscInt turb_spanstats_collect_interval; 127 PetscInt turb_spanstats_viewer_interval; 128 PetscViewer turb_spanstats_viewer; 129 PetscViewerFormat turb_spanstats_viewer_format; 130 // Wall forces 131 struct { 132 PetscInt num_wall; 133 PetscInt *walls; 134 PetscViewer viewer; 135 PetscViewerFormat viewer_format; 136 PetscBool header_written; 137 } wall_forces; 138 // Subgrid Stress Model 139 SGSModelType sgs_model_type; 140 PetscBool sgs_train_enable; 141 // Differential Filtering 142 PetscBool diff_filter_monitor; 143 MeshTransformType mesh_transform_type; 144 }; 145 146 // libCEED data struct 147 struct CeedData_private { 148 CeedVector x_coord, q_data; 149 CeedBasis basis_x, basis_q; 150 CeedElemRestriction elem_restr_x, elem_restr_q, elem_restr_qd_i; 151 OperatorApplyContext op_ics_ctx; 152 }; 153 154 typedef struct { 155 DM dm; 156 PetscSF sf; // For communicating child data to parents 157 OperatorApplyContext op_stats_collect_ctx, op_proj_rhs_ctx; 158 PetscInt num_comp_stats; 159 Vec Child_Stats_loc, Parent_Stats_loc; 160 KSP ksp; // For the L^2 projection solve 161 CeedScalar span_width; // spanwise width of the child domain 162 PetscBool do_mms_test; 163 OperatorApplyContext mms_error_ctx; 164 CeedContextFieldLabel solution_time_label, previous_time_label; 165 } SpanStatsData; 166 167 typedef struct { 168 DM dm; 169 PetscInt num_comp; 170 OperatorApplyContext l2_rhs_ctx; 171 KSP ksp; 172 } *NodalProjectionData; 173 174 typedef PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 175 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 176 typedef struct { 177 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 178 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 179 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 180 CeedVector sgs_nodal_ceed, grad_velo_ceed; 181 SgsDDNodalStressEval sgs_nodal_eval; 182 SgsDDNodalStressInference sgs_nodal_inference; 183 void *sgs_nodal_inference_ctx; 184 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 185 } *SgsDDData; 186 187 typedef struct { 188 DM dm_dd_training; 189 PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 190 PetscScalar filter_widths[16]; 191 OperatorApplyContext op_training_data_calc_ctx; 192 NodalProjectionData filtered_grad_velo_proj; 193 size_t training_data_array_dims[2]; 194 PetscBool overwrite_training_data; 195 } *SGS_DD_TrainingData; 196 197 typedef struct { 198 DM dm_filter; 199 PetscInt num_filtered_fields; 200 CeedInt *num_field_components; 201 PetscInt field_prim_state, field_velo_prod; 202 OperatorApplyContext op_rhs_ctx; 203 KSP ksp; 204 PetscObjectState X_loc_state; 205 PetscBool do_mms_test; 206 CeedContextFieldLabel filter_width_scaling_label; 207 } *DiffFilterData; 208 209 typedef struct { 210 void *client; 211 char rank_id_name[16]; 212 PetscInt collocated_database_num_ranks; 213 } *SmartSimData; 214 215 // PETSc user data 216 struct User_private { 217 MPI_Comm comm; 218 DM dm; 219 DM dm_viz; 220 Mat interp_viz; 221 Ceed ceed; 222 Units units; 223 Vec Q_loc, Q_dot_loc; 224 Physics phys; 225 AppCtx app_ctx; 226 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 227 CeedOperator op_ifunction; 228 Mat mat_ijacobian; 229 KSP mass_ksp; 230 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 231 CeedScalar time_bc_set; 232 SpanStatsData spanstats; 233 NodalProjectionData grad_velo_proj; 234 SgsDDData sgs_dd_data; 235 DiffFilterData diff_filter; 236 SmartSimData smartsim; 237 SGS_DD_TrainingData sgs_dd_train; 238 }; 239 240 // Units 241 struct Units_private { 242 // fundamental units 243 PetscScalar meter; 244 PetscScalar kilogram; 245 PetscScalar second; 246 PetscScalar Kelvin; 247 // derived units 248 PetscScalar Pascal; 249 PetscScalar J_per_kg_K; 250 PetscScalar m_per_squared_s; 251 PetscScalar W_per_m_K; 252 PetscScalar Joule; 253 }; 254 255 // Boundary conditions 256 struct SimpleBC_private { 257 PetscInt num_wall, // Number of faces with wall BCs 258 wall_comps[5], // An array of constrained component numbers 259 num_comps, 260 num_symmetry[3], // Number of faces with symmetry BCs 261 num_inflow, num_outflow, num_freestream, num_slip; 262 PetscInt walls[16], symmetries[3][16], inflows[16], outflows[16], freestreams[16], slips[16]; 263 }; 264 265 // Struct that contains all enums and structs used for the physics of all problems 266 struct Physics_private { 267 PetscBool implicit; 268 StateVariable state_var; 269 CeedContextFieldLabel solution_time_label; 270 CeedContextFieldLabel stg_solution_time_label; 271 CeedContextFieldLabel timestep_size_label; 272 CeedContextFieldLabel ics_time_label; 273 }; 274 275 typedef struct { 276 CeedQFunctionUser qfunction; 277 const char *qfunction_loc; 278 CeedQFunctionContext qfunction_context; 279 } ProblemQFunctionSpec; 280 281 // Problem specific data 282 typedef struct ProblemData_private *ProblemData; 283 struct ProblemData_private { 284 CeedInt dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur; 285 CeedScalar dm_scale; 286 ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip, 287 apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian; 288 bool non_zero_time; 289 PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian; 290 PetscErrorCode (*print_info)(User, ProblemData, AppCtx); 291 PetscErrorCode (*create_mass_operator)(User, CeedOperator *); 292 }; 293 294 extern int FreeContextPetsc(void *); 295 296 // ----------------------------------------------------------------------------- 297 // Set up problems 298 // ----------------------------------------------------------------------------- 299 // Set up function for each problem 300 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 301 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 302 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 303 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 304 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 305 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 306 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 307 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 308 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 309 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc); 310 311 // Print function for each problem 312 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx); 313 314 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx); 315 316 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx); 317 318 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx); 319 320 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData problem, AppCtx app_ctx); 321 322 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts); 323 324 // ----------------------------------------------------------------------------- 325 // libCEED functions 326 // ----------------------------------------------------------------------------- 327 // Utility function to create local CEED restriction 328 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field, 329 CeedElemRestriction *elem_restr); 330 331 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 332 CeedElemRestriction *restriction); 333 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 334 CeedElemRestriction *restriction); 335 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 336 PetscInt q_data_size, CeedElemRestriction *restriction); 337 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 338 PetscInt q_data_size, CeedElemRestriction *restriction); 339 340 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis); 341 342 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc); 343 344 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 345 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 346 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size); 347 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 348 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 349 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size); 350 // ----------------------------------------------------------------------------- 351 // Time-stepping functions 352 // ----------------------------------------------------------------------------- 353 // RHS (Explicit time-stepper) function setup 354 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 355 356 // Implicit time-stepper function setup 357 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 358 359 // User provided TS Monitor 360 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 361 362 // TS: Create, setup, and solve 363 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts); 364 365 // Update Boundary Values when time has changed 366 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t); 367 368 // ----------------------------------------------------------------------------- 369 // Setup DM 370 // ----------------------------------------------------------------------------- 371 // Create mesh 372 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 373 374 // Set up DM 375 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys); 376 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 377 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 378 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 379 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 380 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 381 382 // Refine DM for high-order viz 383 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys); 384 385 // ----------------------------------------------------------------------------- 386 // Process command line options 387 // ----------------------------------------------------------------------------- 388 // Register problems to be available on the command line 389 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx); 390 391 // Process general command line options 392 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc); 393 394 // ----------------------------------------------------------------------------- 395 // Miscellaneous utility functions 396 // ----------------------------------------------------------------------------- 397 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 398 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 399 CeedVector *inv_multiplicity); 400 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time); 401 402 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 403 Vec grad_FVM); 404 405 // Compare reference solution values with current test run for CI 406 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 407 408 // Get error for problems with exact solutions 409 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time); 410 411 // Post-processing 412 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time); 413 414 // -- Gather initial Q values in case of continuation of simulation 415 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q); 416 417 // Record boundary values from initial condition 418 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 419 420 // Versioning token for binary checkpoints 421 extern const PetscInt32 FLUIDS_FILE_TOKEN; // for backwards compatibility 422 extern const PetscInt32 FLUIDS_FILE_TOKEN_32; 423 extern const PetscInt32 FLUIDS_FILE_TOKEN_64; 424 425 // Create appropriate mass qfunction based on number of components N 426 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 427 428 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context); 429 430 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 431 FILE **fp); 432 433 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows); 434 435 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]); 436 437 // ----------------------------------------------------------------------------- 438 // Turbulence Statistics Collection Functions 439 // ----------------------------------------------------------------------------- 440 441 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 442 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 443 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data); 444 445 // ----------------------------------------------------------------------------- 446 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 447 // ----------------------------------------------------------------------------- 448 449 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 450 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 451 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc); 452 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input, 453 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 454 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 455 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 456 CeedVector *grid_aniso_vector); 457 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 458 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 459 460 // ----------------------------------------------------------------------------- 461 // Boundary Condition Related Functions 462 // ----------------------------------------------------------------------------- 463 464 // Setup StrongBCs that use QFunctions 465 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc); 466 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 475 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 476 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 477 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 478 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 479 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem); 480 481 // ----------------------------------------------------------------------------- 482 // SGS Data-Driven Training via SmartSim 483 // ----------------------------------------------------------------------------- 484 PetscErrorCode SmartSimSetup(User user); 485 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 486 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem); 487 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 488 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 489 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 490