1 // Copyright (c) 2017-2022, 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 8 #ifndef libceed_fluids_examples_navier_stokes_h 9 #define libceed_fluids_examples_navier_stokes_h 10 11 #include <ceed.h> 12 #include <mat-ceed.h> 13 #include <petscts.h> 14 #include <stdbool.h> 15 16 #include "./include/petsc_ops.h" 17 #include "qfunctions/newtonian_types.h" 18 #include "qfunctions/stabilization_types.h" 19 20 #if PETSC_VERSION_LT(3, 20, 0) 21 #error "PETSc v3.20 or later is required" 22 #endif 23 24 #if PETSC_VERSION_LT(3, 21, 0) 25 #define DMSetCoordinateDisc(a, b, c) DMProjectCoordinates(a, b) 26 #endif 27 28 #define PetscCallCeed(ceed, ...) \ 29 do { \ 30 int ierr = __VA_ARGS__; \ 31 if (ierr != CEED_ERROR_SUCCESS) { \ 32 const char *error_message; \ 33 CeedGetErrorMessage(ceed, &error_message); \ 34 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "%s", error_message); \ 35 } \ 36 } while (0) 37 38 // ----------------------------------------------------------------------------- 39 // Enums 40 // ----------------------------------------------------------------------------- 41 // Translate PetscMemType to CeedMemType 42 static inline CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; } 43 44 // Euler - test cases 45 typedef enum { 46 EULER_TEST_ISENTROPIC_VORTEX = 0, 47 EULER_TEST_1 = 1, 48 EULER_TEST_2 = 2, 49 EULER_TEST_3 = 3, 50 EULER_TEST_4 = 4, 51 EULER_TEST_5 = 5, 52 } EulerTestType; 53 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1", "test_2", "test_3", "test_4", "test_5", 54 "EulerTestType", "EULER_TEST_", NULL}; 55 56 // Advection - Wind types 57 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL}; 58 59 // Advection - Initial Condition Types 60 static const char *const AdvectionICTypes[] = {"sphere", "cylinder", "cosine_hill", "skew", "AdvectionICType", "ADVECTIONIC_", NULL}; 61 62 // Advection - Bubble Continuity Types 63 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "cosine", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL}; 64 65 // Stabilization methods 66 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 67 68 // Stabilization tau constants 69 static const char *const StabilizationTauTypes[] = {"Ctau", "AdvDiff_Shakib", "AdvDiff_Shakib_P", "StabilizationTauType", "STAB_TAU_", NULL}; 70 71 // Test mode type 72 typedef enum { 73 TESTTYPE_NONE = 0, 74 TESTTYPE_SOLVER = 1, 75 TESTTYPE_TURB_SPANSTATS = 2, 76 TESTTYPE_DIFF_FILTER = 3, 77 } TestType; 78 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "diff_filter", "TestType", "TESTTYPE_", NULL}; 79 80 // Subgrid-Stress mode type 81 typedef enum { 82 SGS_MODEL_NONE = 0, 83 SGS_MODEL_DATA_DRIVEN = 1, 84 } SGSModelType; 85 static const char *const SGSModelTypes[] = {"none", "data_driven", "SGSModelType", "SGS_MODEL_", NULL}; 86 87 // Mesh transformation type 88 typedef enum { 89 MESH_TRANSFORM_NONE = 0, 90 MESH_TRANSFORM_PLATEMESH = 1, 91 } MeshTransformType; 92 static const char *const MeshTransformTypes[] = {"none", "platemesh", "MeshTransformType", "MESH_TRANSFORM_", NULL}; 93 94 static const char *const DifferentialFilterDampingFunctions[] = { 95 "none", "van_driest", "mms", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL}; 96 97 // ----------------------------------------------------------------------------- 98 // Log Events 99 // ----------------------------------------------------------------------------- 100 extern PetscLogEvent FLUIDS_CeedOperatorApply; 101 extern PetscLogEvent FLUIDS_CeedOperatorAssemble; 102 extern PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 103 extern PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 104 extern PetscLogEvent FLUIDS_SmartRedis_Init; 105 extern PetscLogEvent FLUIDS_SmartRedis_Meta; 106 extern PetscLogEvent FLUIDS_SmartRedis_Train; 107 extern PetscLogEvent FLUIDS_TrainDataCompute; 108 extern PetscLogEvent FLUIDS_DifferentialFilter; 109 extern PetscLogEvent FLUIDS_VelocityGradientProjection; 110 PetscErrorCode RegisterLogEvents(); 111 112 // ----------------------------------------------------------------------------- 113 // Structs 114 // ----------------------------------------------------------------------------- 115 // Structs declarations 116 typedef struct AppCtx_private *AppCtx; 117 typedef struct CeedData_private *CeedData; 118 typedef struct User_private *User; 119 typedef struct Units_private *Units; 120 typedef struct SimpleBC_private *SimpleBC; 121 typedef struct Physics_private *Physics; 122 123 // Application context from user command line options 124 struct AppCtx_private { 125 // libCEED arguments 126 char ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend 127 PetscInt degree; 128 PetscInt q_extra; 129 // Solver arguments 130 MatType amat_type; 131 // Post-processing arguments 132 PetscInt checkpoint_interval; 133 PetscInt viz_refine; 134 PetscInt cont_steps; 135 PetscReal cont_time; 136 char cont_file[PETSC_MAX_PATH_LEN]; 137 char cont_time_file[PETSC_MAX_PATH_LEN]; 138 char output_dir[PETSC_MAX_PATH_LEN]; 139 PetscBool add_stepnum2bin; 140 PetscBool checkpoint_vtk; 141 // Problem type arguments 142 PetscFunctionList problems; 143 char problem_name[PETSC_MAX_PATH_LEN]; 144 // Test mode arguments 145 TestType test_type; 146 PetscScalar test_tol; 147 char test_file_path[PETSC_MAX_PATH_LEN]; 148 // Turbulent spanwise statistics 149 PetscBool turb_spanstats_enable; 150 PetscInt turb_spanstats_collect_interval; 151 PetscInt turb_spanstats_viewer_interval; 152 PetscViewer turb_spanstats_viewer; 153 PetscViewerFormat turb_spanstats_viewer_format; 154 // Wall forces 155 struct { 156 PetscInt num_wall; 157 PetscInt *walls; 158 PetscViewer viewer; 159 PetscViewerFormat viewer_format; 160 PetscBool header_written; 161 } wall_forces; 162 // Subgrid Stress Model 163 SGSModelType sgs_model_type; 164 PetscBool sgs_train_enable; 165 // Differential Filtering 166 PetscBool diff_filter_monitor; 167 MeshTransformType mesh_transform_type; 168 }; 169 170 // libCEED data struct 171 struct CeedData_private { 172 CeedVector x_coord, q_data; 173 CeedBasis basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur, basis_xc_sur; 174 CeedElemRestriction elem_restr_x, elem_restr_q, elem_restr_qd_i; 175 CeedOperator op_setup_vol; 176 OperatorApplyContext op_ics_ctx; 177 CeedQFunction qf_setup_vol, qf_ics, qf_rhs_vol, qf_ifunction_vol, qf_setup_sur, qf_apply_inflow, qf_apply_inflow_jacobian, qf_apply_outflow, 178 qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian, qf_apply_slip, qf_apply_slip_jacobian; 179 }; 180 181 typedef struct { 182 DM dm; 183 PetscSF sf; // For communicating child data to parents 184 OperatorApplyContext op_stats_collect_ctx, op_proj_rhs_ctx; 185 PetscInt num_comp_stats; 186 Vec Child_Stats_loc, Parent_Stats_loc; 187 KSP ksp; // For the L^2 projection solve 188 CeedScalar span_width; // spanwise width of the child domain 189 PetscBool do_mms_test; 190 OperatorApplyContext mms_error_ctx; 191 CeedContextFieldLabel solution_time_label, previous_time_label; 192 } SpanStatsData; 193 194 typedef struct { 195 DM dm; 196 PetscInt num_comp; 197 OperatorApplyContext l2_rhs_ctx; 198 KSP ksp; 199 } *NodalProjectionData; 200 201 typedef PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 202 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 203 typedef struct { 204 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 205 PetscInt num_comp_sgs; 206 CeedInt num_comp_inputs, num_comp_outputs; 207 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 208 CeedVector sgs_nodal_ceed, grad_velo_ceed; 209 SgsDDNodalStressEval sgs_nodal_eval; 210 SgsDDNodalStressInference sgs_nodal_inference; 211 void *sgs_nodal_inference_ctx; 212 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 213 } *SgsDDData; 214 215 typedef struct { 216 DM dm_dd_training; 217 PetscInt num_comp_dd_inputs, write_data_interval; 218 OperatorApplyContext op_training_data_calc_ctx; 219 NodalProjectionData filtered_grad_velo_proj; 220 size_t training_data_array_dims[2]; 221 PetscBool overwrite_training_data; 222 } *SGS_DD_TrainingData; 223 224 typedef struct { 225 DM dm_filter; 226 PetscInt num_filtered_fields; 227 CeedInt *num_field_components; 228 PetscInt field_prim_state, field_velo_prod; 229 OperatorApplyContext op_rhs_ctx; 230 KSP ksp; 231 PetscBool do_mms_test; 232 } *DiffFilterData; 233 234 typedef struct { 235 void *client; 236 char rank_id_name[16]; 237 PetscInt collocated_database_num_ranks; 238 } *SmartSimData; 239 240 // PETSc user data 241 struct User_private { 242 MPI_Comm comm; 243 DM dm; 244 DM dm_viz; 245 Mat interp_viz; 246 Ceed ceed; 247 Units units; 248 Vec Q_loc, Q_dot_loc; 249 Physics phys; 250 AppCtx app_ctx; 251 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 252 CeedOperator op_rhs_vol, op_ifunction_vol, op_ifunction; 253 Mat mat_ijacobian; 254 KSP mass_ksp; 255 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 256 CeedScalar time_bc_set; 257 SpanStatsData spanstats; 258 NodalProjectionData grad_velo_proj; 259 SgsDDData sgs_dd_data; 260 DiffFilterData diff_filter; 261 SmartSimData smartsim; 262 SGS_DD_TrainingData sgs_dd_train; 263 }; 264 265 // Units 266 struct Units_private { 267 // fundamental units 268 PetscScalar meter; 269 PetscScalar kilogram; 270 PetscScalar second; 271 PetscScalar Kelvin; 272 // derived units 273 PetscScalar Pascal; 274 PetscScalar J_per_kg_K; 275 PetscScalar m_per_squared_s; 276 PetscScalar W_per_m_K; 277 PetscScalar Joule; 278 }; 279 280 // Boundary conditions 281 struct SimpleBC_private { 282 PetscInt num_wall, // Number of faces with wall BCs 283 wall_comps[5], // An array of constrained component numbers 284 num_comps, 285 num_symmetry[3], // Number of faces with symmetry BCs 286 num_inflow, num_outflow, num_freestream, num_slip; 287 PetscInt walls[16], symmetries[3][16], inflows[16], outflows[16], freestreams[16], slips[16]; 288 }; 289 290 // Struct that contains all enums and structs used for the physics of all problems 291 struct Physics_private { 292 PetscBool implicit; 293 StateVariable state_var; 294 CeedContextFieldLabel solution_time_label; 295 CeedContextFieldLabel stg_solution_time_label; 296 CeedContextFieldLabel timestep_size_label; 297 CeedContextFieldLabel ics_time_label; 298 CeedContextFieldLabel ijacobian_time_shift_label; 299 }; 300 301 typedef struct { 302 CeedQFunctionUser qfunction; 303 const char *qfunction_loc; 304 CeedQFunctionContext qfunction_context; 305 } ProblemQFunctionSpec; 306 307 // Problem specific data 308 typedef struct ProblemData_private ProblemData; 309 struct ProblemData_private { 310 CeedInt dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur; 311 CeedScalar dm_scale; 312 ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, 313 apply_freestream, apply_slip, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian; 314 bool non_zero_time; 315 PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian; 316 PetscErrorCode (*print_info)(User, ProblemData *, AppCtx); 317 }; 318 319 extern int FreeContextPetsc(void *); 320 321 // ----------------------------------------------------------------------------- 322 // Set up problems 323 // ----------------------------------------------------------------------------- 324 // Set up function for each problem 325 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 326 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 327 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 328 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 329 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 330 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 331 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 332 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 333 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 334 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 335 336 // Print function for each problem 337 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData *problem, AppCtx app_ctx); 338 339 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData *problem, AppCtx app_ctx); 340 341 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData *problem, AppCtx app_ctx); 342 343 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx); 344 345 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData *problem, AppCtx app_ctx); 346 347 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm); 348 349 // ----------------------------------------------------------------------------- 350 // libCEED functions 351 // ----------------------------------------------------------------------------- 352 // Utility function to create local CEED restriction 353 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field, 354 CeedElemRestriction *elem_restr); 355 356 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 357 CeedElemRestriction *restriction); 358 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 359 CeedElemRestriction *restriction); 360 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 361 PetscInt q_data_size, CeedElemRestriction *restriction); 362 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 363 PetscInt q_data_size, CeedElemRestriction *restriction); 364 365 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis); 366 367 // Utility function to create CEED Composite Operator for the entire domain 368 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol, 369 CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur, 370 CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian); 371 372 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc); 373 374 // ----------------------------------------------------------------------------- 375 // Time-stepping functions 376 // ----------------------------------------------------------------------------- 377 // Create KSP to solve the inverse mass operator for explicit time stepping schemes 378 PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data); 379 380 // RHS (Explicit time-stepper) function setup 381 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 382 383 // Implicit time-stepper function setup 384 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 385 386 // User provided TS Monitor 387 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 388 389 // TS: Create, setup, and solve 390 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts); 391 392 // Update Boundary Values when time has changed 393 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t); 394 395 // ----------------------------------------------------------------------------- 396 // Setup DM 397 // ----------------------------------------------------------------------------- 398 // Create mesh 399 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm); 400 401 // Set up DM 402 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys); 403 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 404 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 405 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 406 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 407 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 408 409 // Refine DM for high-order viz 410 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys); 411 412 // ----------------------------------------------------------------------------- 413 // Process command line options 414 // ----------------------------------------------------------------------------- 415 // Register problems to be available on the command line 416 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx); 417 418 // Process general command line options 419 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc); 420 421 // ----------------------------------------------------------------------------- 422 // Miscellaneous utility functions 423 // ----------------------------------------------------------------------------- 424 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time); 425 426 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 427 Vec grad_FVM); 428 429 // Compare reference solution values with current test run for CI 430 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 431 432 // Get error for problems with exact solutions 433 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time); 434 435 // Post-processing 436 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time); 437 438 // -- Gather initial Q values in case of continuation of simulation 439 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q); 440 441 // Record boundary values from initial condition 442 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 443 444 // Versioning token for binary checkpoints 445 extern const PetscInt32 FLUIDS_FILE_TOKEN; // for backwards compatibility 446 extern const PetscInt32 FLUIDS_FILE_TOKEN_32; 447 extern const PetscInt32 FLUIDS_FILE_TOKEN_64; 448 449 // Create appropriate mass qfunction based on number of components N 450 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 451 452 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context); 453 454 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 455 FILE **fp); 456 457 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows); 458 459 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]); 460 461 PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc); 462 PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed); 463 464 // ----------------------------------------------------------------------------- 465 // Turbulence Statistics Collection Functions 466 // ----------------------------------------------------------------------------- 467 468 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem); 469 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 470 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data); 471 472 // ----------------------------------------------------------------------------- 473 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 474 // ----------------------------------------------------------------------------- 475 476 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem); 477 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 478 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc); 479 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, StateVariable state_var_input, 480 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 481 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 482 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 483 CeedVector *grid_aniso_vector); 484 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso, 485 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 486 487 // ----------------------------------------------------------------------------- 488 // Boundary Condition Related Functions 489 // ----------------------------------------------------------------------------- 490 491 // Setup StrongBCs that use QFunctions 492 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData *problem, SimpleBC bc); 493 494 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 495 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 496 PetscErrorCode SlipBCSetup(ProblemData *problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 497 498 // ----------------------------------------------------------------------------- 499 // Differential Filtering Functions 500 // ----------------------------------------------------------------------------- 501 502 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem); 503 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 504 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 505 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 506 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData *problem); 507 508 // ----------------------------------------------------------------------------- 509 // SGS Data-Driven Training via SmartSim 510 // ----------------------------------------------------------------------------- 511 PetscErrorCode SmartSimSetup(User user); 512 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 513 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem); 514 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 515 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 516 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 517 518 #endif // libceed_fluids_examples_navier_stokes_h 519