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