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