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