xref: /honee/src/differential_filter.c (revision da59998647b348cb577ad1134294158b446d82cc)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 /// @file
4 /// Functions for setting up and performing differential filtering
5 
6 #include "../qfunctions/differential_filter.h"
7 #include <ceed.h>
8 #include <differential_filter.h>
9 #include <navierstokes.h>
10 #include <petscdmplex.h>
11 
12 const char *const DifferentialFilterDampingFunctions[] = {"NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_",
13                                                           NULL};
14 
15 // @brief Create RHS and LHS operators for differential filtering
16 static PetscErrorCode DifferentialFilterCreateOperators(Honee honee, CeedQFunctionContext diff_filter_qfctx, DiffFilterData diff_filter) {
17   Ceed                ceed      = honee->ceed;
18   DM                  dm_filter = diff_filter->dm_filter;
19   CeedInt             q_data_size;
20   PetscInt            num_comp_x, num_comp_q;
21   CeedVector          q_data;
22   CeedElemRestriction elem_restr_qd;
23   PetscInt            height = 0;
24 
25   PetscFunctionBeginUser;
26   PetscCall(DMGetCoordinateNumComps(honee->dm, &num_comp_x));
27   PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
28 
29   PetscCall(QDataGet(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd,
30                      &q_data, &q_data_size));
31 
32   {  // -- Create RHS MatopApplyContext
33     CeedQFunction       qf_rhs;
34     CeedOperator        op_rhs;
35     CeedBasis           basis_q;
36     CeedElemRestriction elem_restr_q;
37     switch (honee->phys->state_var) {
38       case STATEVAR_PRIMITIVE:
39         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_RHS_Prim, DifferentialFilter_RHS_Prim_loc, &qf_rhs));
40         break;
41       case STATEVAR_CONSERVATIVE:
42         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_RHS_Conserv, DifferentialFilter_RHS_Conserv_loc, &qf_rhs));
43         break;
44       case STATEVAR_ENTROPY:
45         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_RHS_Entropy, DifferentialFilter_RHS_Entropy_loc, &qf_rhs));
46         break;
47     }
48     if (diff_filter->do_mms_test) {
49       PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
50       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_MMS_RHS, DifferentialFilter_MMS_RHS_loc, &qf_rhs));
51     }
52     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
53     PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
54 
55     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, diff_filter_qfctx));
56     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
57     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
58     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "x", num_comp_x, CEED_EVAL_INTERP));
59     for (PetscInt i = 0; i < diff_filter->num_filtered_fields; i++) {
60       char field_name[PETSC_MAX_PATH_LEN];
61 
62       PetscCall(PetscSNPrintf(field_name, PETSC_MAX_PATH_LEN, "v%" PetscInt_FMT, i));
63       PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, field_name, diff_filter->num_field_components[i], CEED_EVAL_INTERP));
64     }
65 
66     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op_rhs));
67     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
68     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
69     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
70     for (PetscInt dm_field = 0; dm_field < diff_filter->num_filtered_fields; dm_field++) {
71       char                field_name[PETSC_MAX_PATH_LEN];
72       CeedElemRestriction elem_restr_filter;
73       CeedBasis           basis_filter;
74 
75       PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_filter));
76       PetscCall(DMPlexCeedBasisCreate(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_filter));
77 
78       PetscCall(PetscSNPrintf(field_name, PETSC_MAX_PATH_LEN, "v%" PetscInt_FMT, dm_field));
79       PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, field_name, elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
80 
81       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filter));
82       PetscCallCeed(ceed, CeedBasisDestroy(&basis_filter));
83     }
84 
85     PetscCall(OperatorApplyContextCreate(honee->dm, dm_filter, ceed, op_rhs, NULL, NULL, honee->Q_loc, NULL, &diff_filter->op_rhs_ctx));
86 
87     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
88     PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
89     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
90     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
91   }
92 
93   {  // Setup LHS Operator and KSP for the differential filtering solve
94     CeedOperator        op_lhs;
95     Mat                 mat_lhs;
96     PetscInt            dim, num_comp_grid_aniso;
97     CeedElemRestriction elem_restr_grid_aniso;
98     CeedVector          grid_aniso_ceed;
99 
100     PetscCall(DMGetDimension(honee->dm, &dim));
101     PetscCall(GridAnisotropyTensorCalculateCollocatedVector(ceed, honee, &elem_restr_grid_aniso, &grid_aniso_ceed, &num_comp_grid_aniso));
102 
103     PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_lhs));
104     for (PetscInt i = 0; i < diff_filter->num_filtered_fields; i++) {
105       CeedQFunction       qf_lhs;
106       PetscInt            num_comp_filter = diff_filter->num_field_components[i];
107       CeedOperator        op_lhs_sub;
108       CeedElemRestriction elem_restr_filter;
109       CeedBasis           basis_filter;
110 
111       switch (num_comp_filter) {
112         case 1:
113           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_1, DifferentialFilter_LHS_1_loc, &qf_lhs));
114           break;
115         case 5:
116           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_5, DifferentialFilter_LHS_5_loc, &qf_lhs));
117           break;
118         case 6:
119           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_6, DifferentialFilter_LHS_6_loc, &qf_lhs));
120           break;
121         case 11:
122           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_11, DifferentialFilter_LHS_11_loc, &qf_lhs));
123           break;
124         default:
125           SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_SUP, "Differential filtering not available for (%" PetscInt_FMT ") components",
126                   num_comp_filter);
127       }
128 
129       PetscCallCeed(ceed, CeedQFunctionSetContext(qf_lhs, diff_filter_qfctx));
130       PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_lhs, 0));
131       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "q", num_comp_filter, CEED_EVAL_INTERP));
132       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "Grad_q", num_comp_filter * dim, CEED_EVAL_GRAD));
133       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
134       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "x", num_comp_x, CEED_EVAL_INTERP));
135       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "qdata", q_data_size, CEED_EVAL_NONE));
136       PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_lhs, "v", num_comp_filter, CEED_EVAL_INTERP));
137       PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_lhs, "Grad_v", num_comp_filter * dim, CEED_EVAL_GRAD));
138 
139       {
140         CeedOperatorField op_field;
141         char              field_name[PETSC_MAX_PATH_LEN];
142 
143         PetscCall(PetscSNPrintf(field_name, PETSC_MAX_PATH_LEN, "v%" PetscInt_FMT, i));
144         PetscCallCeed(ceed, CeedOperatorGetFieldByName(diff_filter->op_rhs_ctx->op, field_name, &op_field));
145         PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filter, &basis_filter, NULL));
146       }
147 
148       PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_lhs, NULL, NULL, &op_lhs_sub));
149       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "q", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
150       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "Grad_q", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
151       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "anisotropy tensor", elem_restr_grid_aniso, CEED_BASIS_NONE, grid_aniso_ceed));
152       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
153       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
154       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "v", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
155       PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "Grad_v", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
156 
157       PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_lhs, op_lhs_sub));
158       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filter));
159       PetscCallCeed(ceed, CeedBasisDestroy(&basis_filter));
160       PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_lhs));
161       PetscCallCeed(ceed, CeedOperatorDestroy(&op_lhs_sub));
162     }
163     PetscCallCeed(ceed, CeedVectorDestroy(&grid_aniso_ceed));
164     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grid_aniso));
165 
166     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_lhs, "filter width scaling", &diff_filter->filter_width_scaling_label));
167     PetscCall(MatCreateCeed(dm_filter, dm_filter, op_lhs, NULL, &mat_lhs));
168 
169     PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm_filter), &diff_filter->ksp));
170     PetscCall(KSPSetOptionsPrefix(diff_filter->ksp, "diff_filter_"));
171     {
172       PC pc;
173       PetscCall(KSPGetPC(diff_filter->ksp, &pc));
174       PetscCall(PCSetType(pc, PCJACOBI));
175       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
176       PetscCall(KSPSetType(diff_filter->ksp, KSPCG));
177       PetscCall(KSPSetNormType(diff_filter->ksp, KSP_NORM_NATURAL));
178       PetscCall(KSPSetTolerances(diff_filter->ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
179     }
180     PetscCall(KSPSetFromOptions_WithMatCeed(diff_filter->ksp, mat_lhs));
181 
182     PetscCall(MatDestroy(&mat_lhs));
183     PetscCallCeed(ceed, CeedOperatorDestroy(&op_lhs));
184   }
185 
186   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
187   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
188   PetscFunctionReturn(PETSC_SUCCESS);
189 }
190 
191 // @brief Setup DM, operators, contexts, etc. for performing differential filtering
192 PetscErrorCode DifferentialFilterSetup(Honee honee, DiffFilterData *diff_filter) {
193   MPI_Comm                  comm    = honee->comm;
194   Ceed                      ceed    = honee->ceed;
195   ProblemData               problem = honee->problem_data;
196   DiffFilterData            diff_filter_;
197   NewtonianIdealGasContext  newt_ctx;
198   DifferentialFilterContext diff_filter_ctx;
199   CeedQFunctionContext      diff_filter_qfctx;
200 
201   PetscFunctionBeginUser;
202   PetscCall(PetscNew(&diff_filter_));
203   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_->do_mms_test, NULL));
204 
205   {  // Create DM for filtered quantities
206     PetscSection section;
207 
208     PetscCall(DMClone(honee->dm, &diff_filter_->dm_filter));
209     PetscCall(DMSetMatrixPreallocateSkip(diff_filter_->dm_filter, PETSC_TRUE));
210     PetscCall(PetscObjectSetName((PetscObject)diff_filter_->dm_filter, "Differential Filtering"));
211 
212     diff_filter_->num_filtered_fields = diff_filter_->do_mms_test ? 1 : 2;
213     PetscCall(PetscMalloc1(diff_filter_->num_filtered_fields, &diff_filter_->num_field_components));
214 
215     if (diff_filter_->do_mms_test) {
216       PetscInt field_components;
217       diff_filter_->num_field_components[0] = field_components = 1;
218       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, diff_filter_->num_filtered_fields,
219                                    &field_components, diff_filter_->dm_filter));
220 
221       PetscCall(DMGetLocalSection(diff_filter_->dm_filter, &section));
222       PetscCall(PetscSectionSetFieldName(section, 0, ""));
223       PetscCall(PetscSectionSetComponentName(section, 0, 0, "FilteredPhi"));
224     } else {
225       PetscInt field_components[2];
226       diff_filter_->num_field_components[0] = field_components[0] = DIFF_FILTER_STATE_NUM;
227       diff_filter_->num_field_components[1] = field_components[1] = DIFF_FILTER_VELOCITY_SQUARED_NUM;
228       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, diff_filter_->num_filtered_fields,
229                                    field_components, diff_filter_->dm_filter));
230 
231       diff_filter_->field_prim_state = 0;
232       diff_filter_->field_velo_prod  = 1;
233       PetscCall(DMGetLocalSection(diff_filter_->dm_filter, &section));
234       PetscCall(PetscSectionSetFieldName(section, diff_filter_->field_prim_state, "Filtered Primitive State Variables"));
235       PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_PRESSURE, "FilteredPressure"));
236       PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_X, "FilteredVelocityX"));
237       PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_Y, "FilteredVelocityY"));
238       PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_Z, "FilteredVelocityZ"));
239       PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_TEMPERATURE, "FilteredTemperature"));
240       PetscCall(PetscSectionSetFieldName(section, diff_filter_->field_velo_prod, "Filtered Velocity Products"));
241       PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XX, "FilteredVelocitySquaredXX"));
242       PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_YY, "FilteredVelocitySquaredYY"));
243       PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_ZZ, "FilteredVelocitySquaredZZ"));
244       PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_YZ, "FilteredVelocitySquaredYZ"));
245       PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XZ, "FilteredVelocitySquaredXZ"));
246       PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XY, "FilteredVelocitySquaredXY"));
247     }
248   }
249 
250   PetscCall(PetscNew(&diff_filter_ctx));
251   *diff_filter_ctx = (struct DifferentialFilterContext_){
252       .grid_based_width = false,
253       .width_scaling    = {1, 1, 1},
254       .kernel_scaling   = 0.1,
255       .damping_function = DIFF_FILTER_DAMP_NONE,
256       .friction_length  = 0,
257       .damping_constant = 25,
258   };
259 
260   PetscOptionsBegin(comm, NULL, "Differential Filtering Options", NULL);
261   PetscInt narray = 3;
262   PetscCall(PetscOptionsBool("-diff_filter_grid_based_width", "Use filter width based on the grid size", NULL, diff_filter_ctx->grid_based_width,
263                              (PetscBool *)&diff_filter_ctx->grid_based_width, NULL));
264   PetscCall(PetscOptionsRealArray("-diff_filter_width_scaling", "Anisotropic scaling of filter width tensor", NULL, diff_filter_ctx->width_scaling,
265                                   &narray, NULL));
266   PetscCall(PetscOptionsReal("-diff_filter_kernel_scaling", "Scaling to make differential kernel size \"equivalent\" to other filter kernels", NULL,
267                              diff_filter_ctx->kernel_scaling, &diff_filter_ctx->kernel_scaling, NULL));
268   PetscCall(PetscOptionsEnum("-diff_filter_wall_damping_function", "Damping function to use at the wall", NULL, DifferentialFilterDampingFunctions,
269                              (PetscEnum)diff_filter_ctx->damping_function, (PetscEnum *)&diff_filter_ctx->damping_function, NULL));
270   PetscCall(PetscOptionsReal("-diff_filter_wall_damping_constant", "Contant for the wall-damping function", NULL, diff_filter_ctx->damping_constant,
271                              &diff_filter_ctx->damping_constant, NULL));
272   PetscCall(PetscOptionsReal("-diff_filter_friction_length", "Friction length associated with the flow, \\delta_\\nu. For wall-damping functions",
273                              NULL, diff_filter_ctx->friction_length, &diff_filter_ctx->friction_length, NULL));
274   PetscOptionsEnd();
275 
276   Units units = honee->units;
277   for (int i = 0; i < 3; i++) diff_filter_ctx->width_scaling[i] *= units->meter;
278   diff_filter_ctx->kernel_scaling *= units->meter;
279   diff_filter_ctx->friction_length *= units->meter;
280 
281   // -- Create QFContext
282   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx));
283   diff_filter_ctx->newt_ctx = *newt_ctx;
284   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx));
285 
286   PetscCallCeed(ceed, CeedQFunctionContextCreate(ceed, &diff_filter_qfctx));
287   PetscCallCeed(ceed, CeedQFunctionContextSetData(diff_filter_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*diff_filter_ctx), diff_filter_ctx));
288   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(diff_filter_qfctx, CEED_MEM_HOST, FreeContextPetsc));
289   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(diff_filter_qfctx, "filter width scaling",
290                                                          offsetof(struct DifferentialFilterContext_, width_scaling),
291                                                          sizeof(diff_filter_ctx->width_scaling) / sizeof(diff_filter_ctx->width_scaling[0]),
292                                                          "Filter width scaling"));
293 
294   // -- Setup Operators
295   PetscCall(DifferentialFilterCreateOperators(honee, diff_filter_qfctx, diff_filter_));
296   *diff_filter = diff_filter_;
297 
298   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&diff_filter_qfctx));
299   PetscFunctionReturn(PETSC_SUCCESS);
300 }
301 
302 // @brief Apply differential filter to the solution given by Q
303 PetscErrorCode DifferentialFilterApply(Honee honee, DiffFilterData diff_filter, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution) {
304   PetscObjectState X_loc_state;
305   Vec              RHS;
306 
307   PetscFunctionBeginUser;
308   PetscCall(PetscLogEventBegin(HONEE_DifferentialFilter, Q, Filtered_Solution, 0, 0));
309   PetscCall(DMGetNamedGlobalVector(diff_filter->dm_filter, "RHS", &RHS));
310   PetscCall(UpdateBoundaryValues(honee, diff_filter->op_rhs_ctx->X_loc, solution_time));
311   PetscCall(VecGetState(diff_filter->op_rhs_ctx->X_loc, &X_loc_state));
312   if (X_loc_state != diff_filter->X_loc_state) {
313     PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, RHS, diff_filter->op_rhs_ctx));
314     PetscCall(VecGetState(diff_filter->op_rhs_ctx->X_loc, &X_loc_state));
315     diff_filter->X_loc_state = X_loc_state;
316   }
317   PetscCall(VecViewFromOptions(RHS, NULL, "-diff_filter_rhs_view"));
318 
319   PetscCall(KSPSolve(diff_filter->ksp, RHS, Filtered_Solution));
320   PetscCall(DMRestoreNamedGlobalVector(diff_filter->dm_filter, "RHS", &RHS));
321   PetscCall(PetscLogEventEnd(HONEE_DifferentialFilter, Q, Filtered_Solution, 0, 0));
322   PetscFunctionReturn(PETSC_SUCCESS);
323 }
324 
325 PetscErrorCode TSMonitor_DifferentialFilterSetup(TS ts, PetscViewerAndFormat *ctx) {
326   Honee          honee;
327   DiffFilterData diff_filter;
328 
329   PetscFunctionBeginUser;
330   PetscCall(TSGetApplicationContext(ts, &honee));
331   PetscCall(DifferentialFilterSetup(honee, &diff_filter));
332   ctx->data         = diff_filter;
333   ctx->data_destroy = (PetscCtxDestroyFn *)DifferentialFilterDataDestroy;
334   PetscFunctionReturn(PETSC_SUCCESS);
335 }
336 
337 // @brief TSMonitor for just applying differential filtering to the simulation
338 // This runs every time step and is primarily for testing purposes
339 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
340   Honee          honee;
341   DiffFilterData diff_filter = (DiffFilterData)ctx->data;
342   Vec            Filtered_Field;
343 
344   PetscFunctionBeginUser;
345   PetscCall(TSGetApplicationContext(ts, &honee));
346   PetscCall(DMGetGlobalVector(diff_filter->dm_filter, &Filtered_Field));
347 
348   PetscCall(DifferentialFilterApply(honee, diff_filter, solution_time, Q, Filtered_Field));
349   PetscCall(VecViewFromOptions(Filtered_Field, NULL, "-diff_filter_view"));
350   if (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER) PetscCall(RegressionTest(honee->app_ctx, Filtered_Field));
351 
352   PetscCall(DMRestoreGlobalVector(diff_filter->dm_filter, &Filtered_Field));
353   PetscFunctionReturn(PETSC_SUCCESS);
354 }
355 
356 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData *diff_filter) {
357   DiffFilterData diff_filter_ = *diff_filter;
358 
359   PetscFunctionBeginUser;
360   if (!diff_filter_) PetscFunctionReturn(PETSC_SUCCESS);
361   PetscCall(OperatorApplyContextDestroy(diff_filter_->op_rhs_ctx));
362   PetscCall(DMDestroy(&diff_filter_->dm_filter));
363   PetscCall(KSPDestroy(&diff_filter_->ksp));
364   PetscCall(PetscFree(diff_filter_->num_field_components));
365   PetscCall(PetscFree(diff_filter_));
366   *diff_filter = NULL;
367   PetscFunctionReturn(PETSC_SUCCESS);
368 }
369 
370 PetscErrorCode DifferentialFilterMmsICSetup(Honee honee) {
371   PetscFunctionBeginUser;
372   honee->problem_data->ics = (HoneeQFSpec){
373       .qf_func_ptr = DifferentialFilter_MMS_IC,
374       .qf_loc      = DifferentialFilter_MMS_IC_loc,
375       .qfctx       = honee->problem_data->ics.qfctx,
376   };
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379