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