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