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