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