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 // Differential Filtering 120 PetscBool diff_filter_monitor; 121 MeshTransformType mesh_transform_type; 122 // Divergence of Diffusive Flux Projection 123 DivDiffFluxProjectionMethod divFdiffproj_method; 124 125 PetscInt check_step_interval; 126 }; 127 128 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData; 129 130 struct DivDiffFluxProjectionData_private { 131 PetscInt num_diff_flux_comps; 132 DivDiffFluxProjectionMethod method; 133 NodalProjectionData projection; 134 135 // CeedOperator Objects 136 CeedElemRestriction elem_restr_div_diff_flux; 137 CeedBasis basis_div_diff_flux; 138 CeedEvalMode eval_mode_div_diff_flux; 139 CeedVector div_diff_flux_ceed; 140 141 // Problem specific setup functions 142 PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *); 143 PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *); 144 145 // Only used for direct method: 146 Vec DivDiffFlux_loc; 147 PetscMemType DivDiffFlux_memtype; 148 PetscBool ceed_vec_has_array; 149 150 // Only used for indirect method: 151 OperatorApplyContext calc_div_diff_flux; 152 }; 153 154 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 155 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 156 typedef struct { 157 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 158 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 159 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 160 CeedVector sgs_nodal_ceed, grad_velo_ceed; 161 SgsDDNodalStressEval sgs_nodal_eval; 162 SgsDDNodalStressInference sgs_nodal_inference; 163 void *sgs_nodal_inference_ctx; 164 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 165 } *SgsDDData; 166 167 typedef struct { 168 DM dm_dd_training; 169 PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 170 PetscScalar filter_widths[16]; 171 OperatorApplyContext op_training_data_calc_ctx; 172 NodalProjectionData filtered_grad_velo_proj; 173 size_t training_data_array_dims[2]; 174 PetscBool overwrite_training_data; 175 } *SGS_DD_TrainingData; 176 177 typedef struct { 178 DM dm_filter; 179 PetscInt num_filtered_fields; 180 CeedInt *num_field_components; 181 PetscInt field_prim_state, field_velo_prod; 182 OperatorApplyContext op_rhs_ctx; 183 KSP ksp; 184 PetscObjectState X_loc_state; 185 PetscBool do_mms_test; 186 CeedContextFieldLabel filter_width_scaling_label; 187 } *DiffFilterData; 188 189 typedef struct { 190 void *client; 191 char rank_id_name[16]; 192 PetscInt collocated_database_num_ranks; 193 } *SmartSimData; 194 195 typedef struct _HoneeOps *HoneeOps; 196 struct _HoneeOps {}; 197 198 PetscErrorCode HoneeInit(MPI_Comm comm, Honee *honee); 199 PetscErrorCode HoneeDestroy(Honee *honee); 200 201 // PETSc user data 202 struct Honee_private { 203 PETSCHEADER(struct _HoneeOps); 204 MPI_Comm comm; 205 DM dm; 206 DM dm_viz; 207 Mat interp_viz; 208 Ceed ceed; 209 Units units; 210 Vec Q_loc, Q_dot_loc; 211 Physics phys; 212 AppCtx app_ctx; 213 CeedVector q_ceed, q_dot_ceed, g_ceed, x_ceed; 214 CeedOperator op_ifunction; 215 Mat mat_ijacobian; 216 KSP mass_ksp; 217 OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx; 218 CeedScalar time_bc_set; 219 SgsDDData sgs_dd_data; 220 DiffFilterData diff_filter; 221 SmartSimData smartsim; 222 SGS_DD_TrainingData sgs_dd_train; 223 DivDiffFluxProjectionData diff_flux_proj; 224 225 ProblemData problem_data; 226 227 // old CeedData 228 CeedVector x_coord; 229 CeedBasis basis_x, basis_q; 230 CeedElemRestriction elem_restr_x, elem_restr_q; 231 OperatorApplyContext op_ics_ctx; 232 233 PetscBool set_poststep; 234 time_t start_time; 235 time_t max_wall_time; 236 PetscInt max_wall_time_interval; 237 }; 238 239 // Units 240 struct Units_private { 241 // fundamental units 242 PetscScalar meter; 243 PetscScalar kilogram; 244 PetscScalar second; 245 PetscScalar Kelvin; 246 // derived units 247 PetscScalar Pascal; 248 PetscScalar J_per_kg_K; 249 PetscScalar m_per_squared_s; 250 PetscScalar W_per_m_K; 251 PetscScalar Joule; 252 }; 253 254 // Struct that contains all enums and structs used for the physics of all problems 255 struct Physics_private { 256 PetscBool implicit; 257 StateVariable state_var; 258 CeedContextFieldLabel solution_time_label; 259 CeedContextFieldLabel stg_solution_time_label; 260 CeedContextFieldLabel timestep_size_label; 261 CeedContextFieldLabel ics_time_label; 262 }; 263 264 typedef struct { 265 Honee honee; 266 CeedInt num_comps_jac_data; 267 CeedQFunctionContext qfctx; 268 void *ctx; 269 PetscCtxDestroyFn *DestroyCtx; 270 } *HoneeBCStruct; 271 272 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx); 273 PetscErrorCode HoneeBCDestroy(void **ctx); 274 PetscErrorCode HoneeBCCreateIFunctionQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx, 275 CeedQFunction *qf_ifunc); 276 PetscErrorCode HoneeBCCreateIJacobianQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx, 277 CeedQFunction *qf_ijac); 278 PetscErrorCode HoneeBCAddIFunctionOp(BCDefinition bc_def, DMLabel domain_label, PetscInt label_value, CeedQFunction qf_ifunc, CeedOperator op_ifunc, 279 CeedOperator *sub_op_ifunc); 280 PetscErrorCode HoneeBCAddIJacobianOp(BCDefinition bc_def, CeedOperator sub_op_ifunc, DMLabel domain_label, PetscInt label_value, 281 CeedQFunction qf_ijac, CeedOperator op_ijac); 282 283 typedef struct { 284 CeedQFunctionUser qf_func_ptr; // !< QFunction function pointer 285 const char *qf_loc; // !< Absolute path to QFunction source file 286 CeedQFunctionContext qfctx; // !< QFunctionContext to attach to QFunction 287 } ProblemQFunctionSpec; 288 289 // Problem specific data 290 struct ProblemData_private { 291 CeedInt num_comps_jac_data; 292 ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian; 293 bool compute_exact_solution_error; 294 PetscBool set_bc_from_ics, use_strong_bc_ceed; 295 PetscCount num_bc_defs; 296 BCDefinition *bc_defs; 297 PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx); 298 PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *); 299 }; 300 301 extern int FreeContextPetsc(void *); 302 303 // ----------------------------------------------------------------------------- 304 // Set up problems 305 // ----------------------------------------------------------------------------- 306 // Set up function for each problem 307 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx); 308 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx); 309 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx); 310 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx); 311 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx); 312 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx); 313 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx); 314 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx); 315 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx); 316 317 // Print function for each problem 318 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx); 319 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx); 320 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx); 321 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx); 322 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx); 323 324 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts); 325 326 // ----------------------------------------------------------------------------- 327 // libCEED functions 328 // ----------------------------------------------------------------------------- 329 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem); 330 331 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 332 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 333 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size); 334 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 335 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 336 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size); 337 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size); 338 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord, 339 CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size); 340 PetscErrorCode QDataClearStoredData(); 341 // ----------------------------------------------------------------------------- 342 // Time-stepping functions 343 // ----------------------------------------------------------------------------- 344 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 345 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 346 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 347 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts); 348 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t); 349 350 // ----------------------------------------------------------------------------- 351 // Setup DM 352 // ----------------------------------------------------------------------------- 353 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm); 354 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, Physics phys); 355 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, Physics phys); 356 357 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 358 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 359 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm); 360 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 361 PetscInt num_fields, const PetscInt *field_sizes, DM dm); 362 363 // ----------------------------------------------------------------------------- 364 // Process command line options 365 // ----------------------------------------------------------------------------- 366 PetscErrorCode ProcessCommandLineOptions(Honee honee); 367 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]); 368 369 // ----------------------------------------------------------------------------- 370 // Miscellaneous utility functions 371 // ----------------------------------------------------------------------------- 372 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 373 PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity, 374 CeedVector *inv_multiplicity); 375 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time); 376 377 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 378 Vec grad_FVM); 379 380 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q); 381 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time); 382 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time); 383 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc); 384 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 385 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume); 386 387 // ----------------------------------------------------------------------------- 388 // Turbulence Statistics Collection Functions 389 // ----------------------------------------------------------------------------- 390 PetscErrorCode SpanwiseStatisticsSetup_Turbulence(TS ts, PetscViewerAndFormat *ctx); 391 PetscErrorCode TSMonitor_TurbulenceSpanwiseStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 392 393 // ----------------------------------------------------------------------------- 394 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions 395 // ----------------------------------------------------------------------------- 396 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem); 397 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data); 398 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc); 399 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input, 400 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj); 401 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient); 402 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 403 CeedVector *grid_aniso_vector); 404 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 405 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso); 406 407 // ----------------------------------------------------------------------------- 408 // Boundary Condition Related Functions 409 // ----------------------------------------------------------------------------- 410 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem); 411 PetscErrorCode FreestreamBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, 412 const StatePrimitive *reference); 413 PetscErrorCode OutflowBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, 414 const StatePrimitive *reference); 415 PetscErrorCode SlipBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx); 416 417 // ----------------------------------------------------------------------------- 418 // Differential Filtering Functions 419 // ----------------------------------------------------------------------------- 420 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, ProblemData problem); 421 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter); 422 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 423 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution); 424 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem); 425 426 // ----------------------------------------------------------------------------- 427 // SGS Data-Driven Training via SmartSim 428 // ----------------------------------------------------------------------------- 429 PetscErrorCode SmartSimSetup(Honee honee); 430 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim); 431 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem); 432 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx); 433 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts); 434 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train); 435 436 // ----------------------------------------------------------------------------- 437 // Divergence of Diffusive Flux Projection 438 // ----------------------------------------------------------------------------- 439 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj); 440 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 441 CeedVector *vector, CeedEvalMode *eval_mode); 442 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj); 443 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc); 444 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj); 445 446 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx); 447 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 448 449 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx); 450 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx); 451 452 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx); 453