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
DifferentialFilterCreateOperators(Honee honee,CeedQFunctionContext diff_filter_qfctx,DiffFilterData diff_filter)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, x_coord;
22 CeedElemRestriction elem_restr_qd, elem_restr_x;
23 CeedBasis basis_x;
24 PetscInt height = 0;
25
26 PetscFunctionBeginUser;
27 PetscCall(DMGetCoordinateNumComps(honee->dm, &num_comp_x));
28 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
29 PetscCall(DMPlexCeedCoordinateCreateField(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, &elem_restr_x, &basis_x, &x_coord));
30 PetscCall(QDataGet(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &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", elem_restr_x, basis_x, 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", elem_restr_x, basis_x, 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, CeedElemRestrictionDestroy(&elem_restr_x));
187 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x));
188 PetscCallCeed(ceed, CeedVectorDestroy(&x_coord));
189 PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
190 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
191 PetscFunctionReturn(PETSC_SUCCESS);
192 }
193
194 // @brief Setup DM, operators, contexts, etc. for performing differential filtering
DifferentialFilterSetup(Honee honee,DiffFilterData * diff_filter)195 PetscErrorCode DifferentialFilterSetup(Honee honee, DiffFilterData *diff_filter) {
196 MPI_Comm comm = honee->comm;
197 Ceed ceed = honee->ceed;
198 ProblemData problem = honee->problem_data;
199 DiffFilterData diff_filter_;
200 NewtonianIdealGasContext newt_ctx;
201 DifferentialFilterContext diff_filter_ctx;
202 CeedQFunctionContext diff_filter_qfctx;
203
204 PetscFunctionBeginUser;
205 PetscCall(PetscNew(&diff_filter_));
206 PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_->do_mms_test, NULL));
207
208 { // Create DM for filtered quantities
209 PetscSection section;
210
211 PetscCall(DMClone(honee->dm, &diff_filter_->dm_filter));
212 PetscCall(DMSetMatrixPreallocateSkip(diff_filter_->dm_filter, PETSC_TRUE));
213 PetscCall(PetscObjectSetName((PetscObject)diff_filter_->dm_filter, "Differential Filtering"));
214
215 diff_filter_->num_filtered_fields = diff_filter_->do_mms_test ? 1 : 2;
216 PetscCall(PetscMalloc1(diff_filter_->num_filtered_fields, &diff_filter_->num_field_components));
217
218 if (diff_filter_->do_mms_test) {
219 PetscInt field_components;
220 diff_filter_->num_field_components[0] = field_components = 1;
221 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, diff_filter_->num_filtered_fields,
222 &field_components, diff_filter_->dm_filter));
223
224 PetscCall(DMGetLocalSection(diff_filter_->dm_filter, §ion));
225 PetscCall(PetscSectionSetFieldName(section, 0, ""));
226 PetscCall(PetscSectionSetComponentName(section, 0, 0, "FilteredPhi"));
227 } else {
228 PetscInt field_components[2];
229 diff_filter_->num_field_components[0] = field_components[0] = DIFF_FILTER_STATE_NUM;
230 diff_filter_->num_field_components[1] = field_components[1] = DIFF_FILTER_VELOCITY_SQUARED_NUM;
231 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, diff_filter_->num_filtered_fields,
232 field_components, diff_filter_->dm_filter));
233
234 diff_filter_->field_prim_state = 0;
235 diff_filter_->field_velo_prod = 1;
236 PetscCall(DMGetLocalSection(diff_filter_->dm_filter, §ion));
237 PetscCall(PetscSectionSetFieldName(section, diff_filter_->field_prim_state, "Filtered Primitive State Variables"));
238 PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_PRESSURE, "FilteredPressure"));
239 PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_X, "FilteredVelocityX"));
240 PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_Y, "FilteredVelocityY"));
241 PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_Z, "FilteredVelocityZ"));
242 PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_TEMPERATURE, "FilteredTemperature"));
243 PetscCall(PetscSectionSetFieldName(section, diff_filter_->field_velo_prod, "Filtered Velocity Products"));
244 PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XX, "FilteredVelocitySquaredXX"));
245 PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_YY, "FilteredVelocitySquaredYY"));
246 PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_ZZ, "FilteredVelocitySquaredZZ"));
247 PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_YZ, "FilteredVelocitySquaredYZ"));
248 PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XZ, "FilteredVelocitySquaredXZ"));
249 PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XY, "FilteredVelocitySquaredXY"));
250 }
251 }
252
253 PetscCall(PetscNew(&diff_filter_ctx));
254 *diff_filter_ctx = (struct DifferentialFilterContext_){
255 .grid_based_width = false,
256 .width_scaling = {1, 1, 1},
257 .kernel_scaling = 0.1,
258 .damping_function = DIFF_FILTER_DAMP_NONE,
259 .friction_length = 0,
260 .damping_constant = 25,
261 };
262
263 PetscOptionsBegin(comm, NULL, "Differential Filtering Options", NULL);
264 PetscInt narray = 3;
265 PetscCall(PetscOptionsBool("-diff_filter_grid_based_width", "Use filter width based on the grid size", NULL, diff_filter_ctx->grid_based_width,
266 (PetscBool *)&diff_filter_ctx->grid_based_width, NULL));
267 PetscCall(PetscOptionsRealArray("-diff_filter_width_scaling", "Anisotropic scaling of filter width tensor", NULL, diff_filter_ctx->width_scaling,
268 &narray, NULL));
269 PetscCall(PetscOptionsReal("-diff_filter_kernel_scaling", "Scaling to make differential kernel size \"equivalent\" to other filter kernels", NULL,
270 diff_filter_ctx->kernel_scaling, &diff_filter_ctx->kernel_scaling, NULL));
271 PetscCall(PetscOptionsEnum("-diff_filter_wall_damping_function", "Damping function to use at the wall", NULL, DifferentialFilterDampingFunctions,
272 (PetscEnum)diff_filter_ctx->damping_function, (PetscEnum *)&diff_filter_ctx->damping_function, NULL));
273 PetscCall(PetscOptionsReal("-diff_filter_wall_damping_constant", "Contant for the wall-damping function", NULL, diff_filter_ctx->damping_constant,
274 &diff_filter_ctx->damping_constant, NULL));
275 PetscCall(PetscOptionsReal("-diff_filter_friction_length", "Friction length associated with the flow, \\delta_\\nu. For wall-damping functions",
276 NULL, diff_filter_ctx->friction_length, &diff_filter_ctx->friction_length, NULL));
277 PetscOptionsEnd();
278
279 Units units = honee->units;
280 for (int i = 0; i < 3; i++) diff_filter_ctx->width_scaling[i] *= units->meter;
281 diff_filter_ctx->kernel_scaling *= units->meter;
282 diff_filter_ctx->friction_length *= units->meter;
283
284 // -- Create QFContext
285 PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx));
286 diff_filter_ctx->newt_ctx = *newt_ctx;
287 PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx));
288
289 PetscCallCeed(ceed, CeedQFunctionContextCreate(ceed, &diff_filter_qfctx));
290 PetscCallCeed(ceed, CeedQFunctionContextSetData(diff_filter_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*diff_filter_ctx), diff_filter_ctx));
291 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(diff_filter_qfctx, CEED_MEM_HOST, FreeContextPetsc));
292 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(diff_filter_qfctx, "filter width scaling",
293 offsetof(struct DifferentialFilterContext_, width_scaling),
294 sizeof(diff_filter_ctx->width_scaling) / sizeof(diff_filter_ctx->width_scaling[0]),
295 "Filter width scaling"));
296
297 // -- Setup Operators
298 PetscCall(DifferentialFilterCreateOperators(honee, diff_filter_qfctx, diff_filter_));
299 *diff_filter = diff_filter_;
300
301 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&diff_filter_qfctx));
302 PetscFunctionReturn(PETSC_SUCCESS);
303 }
304
305 // @brief Apply differential filter to the solution given by Q
DifferentialFilterApply(Honee honee,DiffFilterData diff_filter,const PetscReal solution_time,const Vec Q,Vec Filtered_Solution)306 PetscErrorCode DifferentialFilterApply(Honee honee, DiffFilterData diff_filter, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution) {
307 PetscObjectState X_loc_state;
308 Vec RHS;
309
310 PetscFunctionBeginUser;
311 PetscCall(PetscLogEventBegin(HONEE_DifferentialFilter, Q, Filtered_Solution, 0, 0));
312 PetscCall(DMGetNamedGlobalVector(diff_filter->dm_filter, "RHS", &RHS));
313 PetscCall(UpdateBoundaryValues(honee, diff_filter->op_rhs_ctx->X_loc, solution_time));
314 PetscCall(VecGetState(diff_filter->op_rhs_ctx->X_loc, &X_loc_state));
315 if (X_loc_state != diff_filter->X_loc_state) {
316 PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, RHS, diff_filter->op_rhs_ctx));
317 PetscCall(VecGetState(diff_filter->op_rhs_ctx->X_loc, &X_loc_state));
318 diff_filter->X_loc_state = X_loc_state;
319 }
320 PetscCall(VecViewFromOptions(RHS, NULL, "-diff_filter_rhs_view"));
321
322 PetscCall(KSPSolve(diff_filter->ksp, RHS, Filtered_Solution));
323 PetscCall(DMRestoreNamedGlobalVector(diff_filter->dm_filter, "RHS", &RHS));
324 PetscCall(PetscLogEventEnd(HONEE_DifferentialFilter, Q, Filtered_Solution, 0, 0));
325 PetscFunctionReturn(PETSC_SUCCESS);
326 }
327
TSMonitor_DifferentialFilterSetup(TS ts,PetscViewerAndFormat * ctx)328 PetscErrorCode TSMonitor_DifferentialFilterSetup(TS ts, PetscViewerAndFormat *ctx) {
329 Honee honee;
330 DiffFilterData diff_filter;
331
332 PetscFunctionBeginUser;
333 PetscCall(TSGetApplicationContext(ts, &honee));
334 PetscCall(DifferentialFilterSetup(honee, &diff_filter));
335 ctx->data = diff_filter;
336 ctx->data_destroy = (PetscCtxDestroyFn *)DifferentialFilterDataDestroy;
337 PetscFunctionReturn(PETSC_SUCCESS);
338 }
339
340 // @brief TSMonitor for just applying differential filtering to the simulation
341 // This runs every time step and is primarily for testing purposes
TSMonitor_DifferentialFilter(TS ts,PetscInt steps,PetscReal solution_time,Vec Q,PetscViewerAndFormat * ctx)342 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
343 Honee honee;
344 DiffFilterData diff_filter = (DiffFilterData)ctx->data;
345 Vec Filtered_Field;
346
347 PetscFunctionBeginUser;
348 PetscCall(TSGetApplicationContext(ts, &honee));
349 PetscCall(DMGetGlobalVector(diff_filter->dm_filter, &Filtered_Field));
350
351 PetscCall(DifferentialFilterApply(honee, diff_filter, solution_time, Q, Filtered_Field));
352 PetscCall(VecViewFromOptions(Filtered_Field, NULL, "-diff_filter_view"));
353 if (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER) PetscCall(RegressionTest(honee->app_ctx, Filtered_Field));
354
355 PetscCall(DMRestoreGlobalVector(diff_filter->dm_filter, &Filtered_Field));
356 PetscFunctionReturn(PETSC_SUCCESS);
357 }
358
DifferentialFilterDataDestroy(DiffFilterData * diff_filter)359 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData *diff_filter) {
360 DiffFilterData diff_filter_ = *diff_filter;
361
362 PetscFunctionBeginUser;
363 if (!diff_filter_) PetscFunctionReturn(PETSC_SUCCESS);
364 PetscCall(OperatorApplyContextDestroy(diff_filter_->op_rhs_ctx));
365 PetscCall(DMDestroy(&diff_filter_->dm_filter));
366 PetscCall(KSPDestroy(&diff_filter_->ksp));
367 PetscCall(PetscFree(diff_filter_->num_field_components));
368 PetscCall(PetscFree(diff_filter_));
369 *diff_filter = NULL;
370 PetscFunctionReturn(PETSC_SUCCESS);
371 }
372
DifferentialFilterMmsICSetup(Honee honee)373 PetscErrorCode DifferentialFilterMmsICSetup(Honee honee) {
374 PetscFunctionBeginUser;
375 honee->problem_data->ics = (HoneeQFSpec){
376 .qf_func_ptr = DifferentialFilter_MMS_IC,
377 .qf_loc = DifferentialFilter_MMS_IC_loc,
378 .qfctx = honee->problem_data->ics.qfctx,
379 };
380 PetscFunctionReturn(PETSC_SUCCESS);
381 }
382