1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 #pragma once 4 5 #include <ceed.h> 6 #include <bc_definition.h> 7 #include <dm-utils.h> 8 #include <honee-file.h> 9 #include <honee.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 #include <time.h> 16 17 #include <nodal_projection.h> 18 19 #include <petsc_ops.h> 20 #include "../qfunctions/newtonian_types.h" 21 22 #if PETSC_VERSION_LT(3, 23, 0) 23 #error "PETSc v3.23 or later is required" 24 #endif 25 26 // ----------------------------------------------------------------------------- 27 // Enums 28 // ----------------------------------------------------------------------------- 29 30 // Euler - test cases 31 typedef enum { 32 EULER_TEST_ISENTROPIC_VORTEX = 0, 33 EULER_TEST_1 = 1, 34 EULER_TEST_2 = 2, 35 EULER_TEST_3 = 3, 36 EULER_TEST_4 = 4, 37 EULER_TEST_5 = 5, 38 } EulerTestType; 39 static const char *const EulerTestTypes[] = {"ISENTROPIC_VORTEX", "1", "2", "3", "4", "5", "EulerTestType", "EULER_TEST_", NULL}; 40 41 // Test mode type 42 typedef enum { 43 TESTTYPE_NONE = 0, 44 TESTTYPE_SOLVER = 1, 45 TESTTYPE_TURB_SPANSTATS = 2, 46 TESTTYPE_DIFF_FILTER = 3, 47 } TestType; 48 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "TestType", "TESTTYPE_", NULL}; 49 50 // Subgrid-Stress mode type 51 typedef enum { 52 SGS_MODEL_NONE = 0, 53 SGS_MODEL_DATA_DRIVEN = 1, 54 } SGSModelType; 55 static const char *const SGSModelTypes[] = {"NONE", "DATA_DRIVEN", "SGSModelType", "SGS_MODEL_", NULL}; 56 57 // Subgrid-Stress mode type 58 typedef enum { 59 SGS_MODEL_DD_FUSED = 0, 60 SGS_MODEL_DD_SEQENTIAL_CEED = 1, 61 SGS_MODEL_DD_SEQENTIAL_TORCH = 2, 62 } SGSModelDDImplementation; 63 static const char *const SGSModelDDImplementations[] = {"FUSED", "SEQUENTIAL_CEED", "SEQUENTIAL_TORCH", "SGSModelDDImplementation", "SGS_MODEL_DD_", 64 NULL}; 65 66 // Mesh transformation type 67 typedef enum { 68 MESH_TRANSFORM_NONE = 0, 69 MESH_TRANSFORM_PLATEMESH = 1, 70 } MeshTransformType; 71 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL}; 72 73 // ----------------------------------------------------------------------------- 74 // Structs 75 // ----------------------------------------------------------------------------- 76 // Structs declarations 77 typedef struct AppCtx_private *AppCtx; 78 typedef struct Units_private *Units; 79 typedef struct SimpleBC_private *SimpleBC; 80 typedef struct Physics_private *Physics; 81 typedef struct ProblemData_private *ProblemData; 82 83 // Application context from user command line options 84 struct AppCtx_private { 85 // libCEED arguments 86 char ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend 87 PetscInt degree; 88 PetscInt q_extra; 89 // Solver arguments 90 MatType amat_type; 91 // Post-processing arguments 92 PetscInt checkpoint_interval; 93 PetscInt viz_refine; 94 PetscBool use_continue_file; 95 PetscInt cont_steps; 96 PetscReal cont_time; 97 char cont_file[PETSC_MAX_PATH_LEN]; 98 char output_dir[PETSC_MAX_PATH_LEN]; 99 PetscBool add_stepnum2bin; 100 PetscBool checkpoint_vtk; 101 // Problem type arguments 102 PetscFunctionList problems; 103 char problem_name[PETSC_MAX_PATH_LEN]; 104 // Test mode arguments 105 TestType test_type; 106 PetscScalar test_tol; 107 char test_file_path[PETSC_MAX_PATH_LEN]; 108 // Wall forces 109 struct { 110 PetscInt num_wall; 111 PetscInt *walls; 112 PetscViewer viewer; 113 PetscViewerFormat viewer_format; 114 PetscBool header_written; 115 } wall_forces; 116 // Subgrid Stress Model 117 SGSModelType sgs_model_type; 118 PetscBool sgs_train_enable; 119 120 MeshTransformType mesh_transform_type; 121 // Divergence of Diffusive Flux Projection 122 DivDiffFluxProjectionMethod divFdiffproj_method; 123 124 PetscInt check_step_interval; 125 }; 126 127 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData; 128 129 struct DivDiffFluxProjectionData_private { 130 PetscInt num_diff_flux_comps; 131 DivDiffFluxProjectionMethod method; 132 NodalProjectionData projection; 133 134 // CeedOperator Objects 135 CeedElemRestriction elem_restr_div_diff_flux; 136 CeedBasis basis_div_diff_flux; 137 CeedEvalMode eval_mode_div_diff_flux; 138 CeedVector div_diff_flux_ceed; 139 140 // Problem specific setup functions 141 PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *); 142 PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *); 143 144 // Only used for direct method: 145 Vec DivDiffFlux_loc; 146 PetscMemType DivDiffFlux_memtype; 147 PetscBool ceed_vec_has_array; 148 149 // Only used for indirect method: 150 OperatorApplyContext calc_div_diff_flux; 151 }; 152 153 typedef struct { 154 void *client; 155 char rank_id_name[16]; 156 PetscInt collocated_database_num_ranks; 157 } *SmartSimData; 158 159 typedef struct _HoneeOps *HoneeOps; 160 struct _HoneeOps {}; 161 162 PetscErrorCode HoneeInit(MPI_Comm comm, Honee *honee); 163 PetscErrorCode HoneeDestroy(Honee *honee); 164 165 // PETSc user data 166 struct Honee_private { 167 PETSCHEADER(struct _HoneeOps); 168 MPI_Comm comm; 169 DM dm; 170 DM dm_viz; 171 Mat interp_viz; 172 Ceed ceed; 173 Units units; 174 Vec Q_loc, Q_dot_loc; 175 Physics phys; 176 AppCtx app_ctx; 177 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 178 CeedOperator op_ifunction; 179 Mat mat_ijacobian; 180 KSP mass_ksp; 181 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 182 CeedScalar time_bc_set; 183 SmartSimData smartsim; 184 DivDiffFluxProjectionData diff_flux_proj; 185 186 ProblemData problem_data; 187 188 // old CeedData 189 CeedVector x_coord; 190 CeedBasis basis_x, basis_q; 191 CeedElemRestriction elem_restr_x, elem_restr_q; 192 OperatorApplyContext op_ics_ctx; 193 194 PetscBool set_poststep; 195 time_t start_time; 196 time_t max_wall_time; 197 PetscInt max_wall_time_interval; 198 }; 199 200 // Units 201 struct Units_private { 202 // fundamental units 203 PetscScalar meter; 204 PetscScalar kilogram; 205 PetscScalar second; 206 PetscScalar Kelvin; 207 // derived units 208 PetscScalar Pascal; 209 PetscScalar J_per_kg_K; 210 PetscScalar m_per_squared_s; 211 PetscScalar W_per_m_K; 212 PetscScalar Joule; 213 }; 214 215 // Struct that contains all enums and structs used for the physics of all problems 216 struct Physics_private { 217 PetscBool implicit; 218 StateVariable state_var; 219 CeedContextFieldLabel solution_time_label; 220 CeedContextFieldLabel stg_solution_time_label; 221 CeedContextFieldLabel timestep_size_label; 222 CeedContextFieldLabel ics_time_label; 223 }; 224 225 typedef struct { 226 Honee honee; 227 CeedInt num_comps_jac_data; 228 CeedQFunctionContext qfctx; 229 void *ctx; 230 PetscCtxDestroyFn *DestroyCtx; 231 } *HoneeBCStruct; 232 233 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx); 234 PetscErrorCode HoneeBCDestroy(void **ctx); 235 PetscErrorCode HoneeBCCreateIFunctionQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx, 236 CeedQFunction *qf_ifunc); 237 PetscErrorCode HoneeBCCreateIJacobianQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx, 238 CeedQFunction *qf_ijac); 239 PetscErrorCode HoneeBCAddIFunctionOp(BCDefinition bc_def, DMLabel domain_label, PetscInt label_value, CeedQFunction qf_ifunc, CeedOperator op_ifunc, 240 CeedOperator *sub_op_ifunc); 241 PetscErrorCode HoneeBCAddIJacobianOp(BCDefinition bc_def, CeedOperator sub_op_ifunc, DMLabel domain_label, PetscInt label_value, 242 CeedQFunction qf_ijac, CeedOperator op_ijac); 243 244 typedef struct { 245 CeedQFunctionUser qf_func_ptr; // !< QFunction function pointer 246 const char *qf_loc; // !< Absolute path to QFunction source file 247 CeedQFunctionContext qfctx; // !< QFunctionContext to attach to QFunction 248 } ProblemQFunctionSpec; 249 250 // Problem specific data 251 struct ProblemData_private { 252 CeedInt num_comps_jac_data; 253 ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian; 254 bool compute_exact_solution_error; 255 PetscBool set_bc_from_ics, use_strong_bc_ceed; 256 PetscCount num_bc_defs; 257 BCDefinition *bc_defs; 258 PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx); 259 PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *); 260 }; 261 262 extern int FreeContextPetsc(void *); 263 264 // ----------------------------------------------------------------------------- 265 // Set up problems 266 // ----------------------------------------------------------------------------- 267 // Set up function for each problem 268 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx); 269 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx); 270 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx); 271 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx); 272 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx); 273 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx); 274 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx); 275 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx); 276 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx); 277 278 // Print function for each problem 279 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx); 280 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx); 281 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx); 282 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx); 283 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx); 284 285 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts); 286 287 // ----------------------------------------------------------------------------- 288 // libCEED functions 289 // ----------------------------------------------------------------------------- 290 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem); 291 292 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 293 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 294 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size); 295 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 296 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 297 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size); 298 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size); 299 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord, 300 CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 301 PetscErrorCode QDataClearStoredData(); 302 // ----------------------------------------------------------------------------- 303 // Time-stepping functions 304 // ----------------------------------------------------------------------------- 305 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 306 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 307 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 308 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts); 309 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t); 310 311 // ----------------------------------------------------------------------------- 312 // Setup DM 313 // ----------------------------------------------------------------------------- 314 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 315 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, Physics phys); 316 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, Physics phys); 317 318 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 319 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 320 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 321 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 322 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 323 324 // ----------------------------------------------------------------------------- 325 // Process command line options 326 // ----------------------------------------------------------------------------- 327 PetscErrorCode ProcessCommandLineOptions(Honee honee); 328 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]); 329 330 // ----------------------------------------------------------------------------- 331 // Miscellaneous utility functions 332 // ----------------------------------------------------------------------------- 333 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 334 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 335 CeedVector *inv_multiplicity); 336 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time); 337 338 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 339 Vec grad_FVM); 340 341 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 342 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time); 343 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time); 344 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 345 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 346 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume); 347 348 // ----------------------------------------------------------------------------- 349 // Turbulence Statistics Collection Functions 350 // ----------------------------------------------------------------------------- 351 PetscErrorCode SpanwiseStatisticsSetup_Turbulence(TS ts, PetscViewerAndFormat *ctx); 352 PetscErrorCode TSMonitor_TurbulenceSpanwiseStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 353 354 // ----------------------------------------------------------------------------- 355 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 356 // ----------------------------------------------------------------------------- 357 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem); 358 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc); 359 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input, 360 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 361 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 362 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 363 CeedVector *grid_aniso_vector); 364 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 365 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 366 367 // ----------------------------------------------------------------------------- 368 // Boundary Condition Related Functions 369 // ----------------------------------------------------------------------------- 370 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem); 371 PetscErrorCode FreestreamBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, 372 const StatePrimitive *reference); 373 PetscErrorCode OutflowBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, 374 const StatePrimitive *reference); 375 PetscErrorCode SlipBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 376 377 // ----------------------------------------------------------------------------- 378 // SGS Data-Driven Training via SmartSim 379 // ----------------------------------------------------------------------------- 380 PetscErrorCode SmartSimSetup(Honee honee); 381 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 382 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem); 383 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 384 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 385 386 // ----------------------------------------------------------------------------- 387 // Divergence of Diffusive Flux Projection 388 // ----------------------------------------------------------------------------- 389 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj); 390 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 391 CeedVector *vector, CeedEvalMode *eval_mode); 392 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj); 393 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc); 394 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj); 395 396 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx); 397 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 398 399 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx); 400 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 401 402 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx); 403