1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
388b07121SJames Wright /// @file
488b07121SJames Wright /// Functions for setting up and performing differential filtering
588b07121SJames Wright
68d78d7c8SJames Wright #include "../qfunctions/differential_filter.h"
787edb941SJames Wright #include <ceed.h>
88d78d7c8SJames Wright #include <differential_filter.h>
9149fb536SJames Wright #include <navierstokes.h>
108d78d7c8SJames Wright #include <petscdmplex.h>
1188b07121SJames Wright
12e5a8cae0SJames Wright const char *const DifferentialFilterDampingFunctions[] = {"NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_",
13e5a8cae0SJames Wright NULL};
14e5a8cae0SJames Wright
1588b07121SJames Wright // @brief Create RHS and LHS operators for differential filtering
DifferentialFilterCreateOperators(Honee honee,CeedQFunctionContext diff_filter_qfctx,DiffFilterData diff_filter)16cb8a476cSJames Wright static PetscErrorCode DifferentialFilterCreateOperators(Honee honee, CeedQFunctionContext diff_filter_qfctx, DiffFilterData diff_filter) {
17cb8a476cSJames Wright Ceed ceed = honee->ceed;
1888b07121SJames Wright DM dm_filter = diff_filter->dm_filter;
19cf8f12c8SJames Wright CeedInt q_data_size;
20cf8f12c8SJames Wright PetscInt num_comp_x, num_comp_q;
21*4fe35dceSJames Wright CeedVector q_data, x_coord;
22*4fe35dceSJames Wright CeedElemRestriction elem_restr_qd, elem_restr_x;
23*4fe35dceSJames Wright CeedBasis basis_x;
24e3db12f8SJames Wright PetscInt height = 0;
2588b07121SJames Wright
2688b07121SJames Wright PetscFunctionBeginUser;
27cf8f12c8SJames Wright PetscCall(DMGetCoordinateNumComps(honee->dm, &num_comp_x));
28cf8f12c8SJames Wright PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
29*4fe35dceSJames Wright PetscCall(DMPlexCeedCoordinateCreateField(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, &elem_restr_x, &basis_x, &x_coord));
309018c49aSJames Wright PetscCall(QDataGet(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
3188b07121SJames Wright
3288b07121SJames Wright { // -- Create RHS MatopApplyContext
3388b07121SJames Wright CeedQFunction qf_rhs;
3488b07121SJames Wright CeedOperator op_rhs;
35cf8f12c8SJames Wright CeedBasis basis_q;
36cf8f12c8SJames Wright CeedElemRestriction elem_restr_q;
370c373b74SJames Wright switch (honee->phys->state_var) {
3888b07121SJames Wright case STATEVAR_PRIMITIVE:
39b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_RHS_Prim, DifferentialFilter_RHS_Prim_loc, &qf_rhs));
4088b07121SJames Wright break;
4188b07121SJames Wright case STATEVAR_CONSERVATIVE:
42b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_RHS_Conserv, DifferentialFilter_RHS_Conserv_loc, &qf_rhs));
4388b07121SJames Wright break;
449b103f75SJames Wright case STATEVAR_ENTROPY:
459b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_RHS_Entropy, DifferentialFilter_RHS_Entropy_loc, &qf_rhs));
469b103f75SJames Wright break;
4788b07121SJames Wright }
4825c92e8fSJames Wright if (diff_filter->do_mms_test) {
49b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
50b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_MMS_RHS, DifferentialFilter_MMS_RHS_loc, &qf_rhs));
5125c92e8fSJames Wright }
52cf8f12c8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
53cf8f12c8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
5425c92e8fSJames Wright
55b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, diff_filter_qfctx));
56b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
57be29160dSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
58b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "x", num_comp_x, CEED_EVAL_INTERP));
59c5e50263SJames Wright for (PetscInt i = 0; i < diff_filter->num_filtered_fields; i++) {
60c5e50263SJames Wright char field_name[PETSC_MAX_PATH_LEN];
61be29160dSJames Wright
62c5e50263SJames Wright PetscCall(PetscSNPrintf(field_name, PETSC_MAX_PATH_LEN, "v%" PetscInt_FMT, i));
63b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, field_name, diff_filter->num_field_components[i], CEED_EVAL_INTERP));
64c5e50263SJames Wright }
6588b07121SJames Wright
66b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op_rhs));
67cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
68be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
69*4fe35dceSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "x", elem_restr_x, basis_x, x_coord));
7015c18037SJames Wright for (PetscInt dm_field = 0; dm_field < diff_filter->num_filtered_fields; dm_field++) {
71c5e50263SJames Wright char field_name[PETSC_MAX_PATH_LEN];
72c5e50263SJames Wright CeedElemRestriction elem_restr_filter;
73c5e50263SJames Wright CeedBasis basis_filter;
74be29160dSJames Wright
75e3db12f8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_filter));
76e3db12f8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, dm_filter, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_filter));
77c5e50263SJames Wright
7815c18037SJames Wright PetscCall(PetscSNPrintf(field_name, PETSC_MAX_PATH_LEN, "v%" PetscInt_FMT, dm_field));
79b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, field_name, elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
8087edb941SJames Wright
8187edb941SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filter));
8287edb941SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_filter));
83c5e50263SJames Wright }
8488b07121SJames Wright
850c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, dm_filter, ceed, op_rhs, NULL, NULL, honee->Q_loc, NULL, &diff_filter->op_rhs_ctx));
8688b07121SJames Wright
87cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
88cf8f12c8SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
89b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
90b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
9188b07121SJames Wright }
9288b07121SJames Wright
9388b07121SJames Wright { // Setup LHS Operator and KSP for the differential filtering solve
9488b07121SJames Wright CeedOperator op_lhs;
953170c09fSJames Wright Mat mat_lhs;
96defe8520SJames Wright PetscInt dim, num_comp_grid_aniso;
9788b07121SJames Wright CeedElemRestriction elem_restr_grid_aniso;
9888b07121SJames Wright CeedVector grid_aniso_ceed;
9988b07121SJames Wright
1000c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim));
101e3663b90SJames Wright PetscCall(GridAnisotropyTensorCalculateCollocatedVector(ceed, honee, &elem_restr_grid_aniso, &grid_aniso_ceed, &num_comp_grid_aniso));
102c5e50263SJames Wright
103da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_lhs));
104c5e50263SJames Wright for (PetscInt i = 0; i < diff_filter->num_filtered_fields; i++) {
105c5e50263SJames Wright CeedQFunction qf_lhs;
106c5e50263SJames Wright PetscInt num_comp_filter = diff_filter->num_field_components[i];
107c5e50263SJames Wright CeedOperator op_lhs_sub;
108c5e50263SJames Wright CeedElemRestriction elem_restr_filter;
109c5e50263SJames Wright CeedBasis basis_filter;
110c5e50263SJames Wright
11188b07121SJames Wright switch (num_comp_filter) {
11288b07121SJames Wright case 1:
113b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_1, DifferentialFilter_LHS_1_loc, &qf_lhs));
11488b07121SJames Wright break;
11588b07121SJames Wright case 5:
116b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_5, DifferentialFilter_LHS_5_loc, &qf_lhs));
11788b07121SJames Wright break;
118c5e50263SJames Wright case 6:
119b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_6, DifferentialFilter_LHS_6_loc, &qf_lhs));
120c5e50263SJames Wright break;
12188b07121SJames Wright case 11:
122b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DifferentialFilter_LHS_11, DifferentialFilter_LHS_11_loc, &qf_lhs));
12388b07121SJames Wright break;
12488b07121SJames Wright default:
1250c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_SUP, "Differential filtering not available for (%" PetscInt_FMT ") components",
126defe8520SJames Wright num_comp_filter);
12788b07121SJames Wright }
12888b07121SJames Wright
129b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_lhs, diff_filter_qfctx));
1303170c09fSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_lhs, 0));
131b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "q", num_comp_filter, CEED_EVAL_INTERP));
132b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "Grad_q", num_comp_filter * dim, CEED_EVAL_GRAD));
133b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
134b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "x", num_comp_x, CEED_EVAL_INTERP));
135be29160dSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_lhs, "qdata", q_data_size, CEED_EVAL_NONE));
136b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_lhs, "v", num_comp_filter, CEED_EVAL_INTERP));
137b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_lhs, "Grad_v", num_comp_filter * dim, CEED_EVAL_GRAD));
13888b07121SJames Wright
139c5e50263SJames Wright {
140c5e50263SJames Wright CeedOperatorField op_field;
141c5e50263SJames Wright char field_name[PETSC_MAX_PATH_LEN];
142be29160dSJames Wright
143c5e50263SJames Wright PetscCall(PetscSNPrintf(field_name, PETSC_MAX_PATH_LEN, "v%" PetscInt_FMT, i));
144b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(diff_filter->op_rhs_ctx->op, field_name, &op_field));
14501e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filter, &basis_filter, NULL));
146c5e50263SJames Wright }
14788b07121SJames Wright
148b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_lhs, NULL, NULL, &op_lhs_sub));
149b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "q", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
150b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "Grad_q", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
15158e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "anisotropy tensor", elem_restr_grid_aniso, CEED_BASIS_NONE, grid_aniso_ceed));
152*4fe35dceSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "x", elem_restr_x, basis_x, x_coord));
153be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
154b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "v", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
155b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_lhs_sub, "Grad_v", elem_restr_filter, basis_filter, CEED_VECTOR_ACTIVE));
156c5e50263SJames Wright
157da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_lhs, op_lhs_sub));
158fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filter));
159fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_filter));
160b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_lhs));
161b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_lhs_sub));
162c5e50263SJames Wright }
16387edb941SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&grid_aniso_ceed));
16487edb941SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grid_aniso));
16587edb941SJames Wright
1666ea7c1aeSJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_lhs, "filter width scaling", &diff_filter->filter_width_scaling_label));
167000d2032SJeremy L Thompson PetscCall(MatCreateCeed(dm_filter, dm_filter, op_lhs, NULL, &mat_lhs));
16888b07121SJames Wright
16988b07121SJames Wright PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm_filter), &diff_filter->ksp));
17088b07121SJames Wright PetscCall(KSPSetOptionsPrefix(diff_filter->ksp, "diff_filter_"));
17188b07121SJames Wright {
17288b07121SJames Wright PC pc;
17388b07121SJames Wright PetscCall(KSPGetPC(diff_filter->ksp, &pc));
17488b07121SJames Wright PetscCall(PCSetType(pc, PCJACOBI));
17588b07121SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
17688b07121SJames Wright PetscCall(KSPSetType(diff_filter->ksp, KSPCG));
17788b07121SJames Wright PetscCall(KSPSetNormType(diff_filter->ksp, KSP_NORM_NATURAL));
17888b07121SJames Wright PetscCall(KSPSetTolerances(diff_filter->ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
17988b07121SJames Wright }
1803170c09fSJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(diff_filter->ksp, mat_lhs));
18188b07121SJames Wright
182b8462d76SJames Wright PetscCall(MatDestroy(&mat_lhs));
183b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_lhs));
18488b07121SJames Wright }
185be29160dSJames Wright
186*4fe35dceSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x));
187*4fe35dceSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x));
188*4fe35dceSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_coord));
189be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
190be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
191d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
19288b07121SJames Wright }
19388b07121SJames Wright
19488b07121SJames Wright // @brief Setup DM, operators, contexts, etc. for performing differential filtering
DifferentialFilterSetup(Honee honee,DiffFilterData * diff_filter)195cb8a476cSJames Wright PetscErrorCode DifferentialFilterSetup(Honee honee, DiffFilterData *diff_filter) {
1960c373b74SJames Wright MPI_Comm comm = honee->comm;
197cb8a476cSJames Wright Ceed ceed = honee->ceed;
198cb8a476cSJames Wright ProblemData problem = honee->problem_data;
199cb8a476cSJames Wright DiffFilterData diff_filter_;
200cde3d787SJames Wright NewtonianIdealGasContext newt_ctx;
20188b07121SJames Wright DifferentialFilterContext diff_filter_ctx;
20288b07121SJames Wright CeedQFunctionContext diff_filter_qfctx;
20388b07121SJames Wright
20488b07121SJames Wright PetscFunctionBeginUser;
205cb8a476cSJames Wright PetscCall(PetscNew(&diff_filter_));
206cb8a476cSJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_->do_mms_test, NULL));
20788b07121SJames Wright
20888b07121SJames Wright { // Create DM for filtered quantities
20988b07121SJames Wright PetscSection section;
21088b07121SJames Wright
211cb8a476cSJames Wright PetscCall(DMClone(honee->dm, &diff_filter_->dm_filter));
212cb8a476cSJames Wright PetscCall(DMSetMatrixPreallocateSkip(diff_filter_->dm_filter, PETSC_TRUE));
213cb8a476cSJames Wright PetscCall(PetscObjectSetName((PetscObject)diff_filter_->dm_filter, "Differential Filtering"));
21488b07121SJames Wright
215cb8a476cSJames Wright diff_filter_->num_filtered_fields = diff_filter_->do_mms_test ? 1 : 2;
216cb8a476cSJames Wright PetscCall(PetscMalloc1(diff_filter_->num_filtered_fields, &diff_filter_->num_field_components));
21788b07121SJames Wright
218cb8a476cSJames Wright if (diff_filter_->do_mms_test) {
2193bd8abf2SJames Wright PetscInt field_components;
220cb8a476cSJames Wright diff_filter_->num_field_components[0] = field_components = 1;
221cb8a476cSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, diff_filter_->num_filtered_fields,
222cb8a476cSJames Wright &field_components, diff_filter_->dm_filter));
22388b07121SJames Wright
224cb8a476cSJames Wright PetscCall(DMGetLocalSection(diff_filter_->dm_filter, §ion));
22588b07121SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, ""));
22625c92e8fSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "FilteredPhi"));
22725c92e8fSJames Wright } else {
228cf171902SJames Wright PetscInt field_components[2];
229cb8a476cSJames Wright diff_filter_->num_field_components[0] = field_components[0] = DIFF_FILTER_STATE_NUM;
230cb8a476cSJames Wright diff_filter_->num_field_components[1] = field_components[1] = DIFF_FILTER_VELOCITY_SQUARED_NUM;
231cb8a476cSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, diff_filter_->num_filtered_fields,
232cb8a476cSJames Wright field_components, diff_filter_->dm_filter));
2331ed803a2SJames Wright
234cb8a476cSJames Wright diff_filter_->field_prim_state = 0;
235cb8a476cSJames Wright diff_filter_->field_velo_prod = 1;
236cb8a476cSJames Wright PetscCall(DMGetLocalSection(diff_filter_->dm_filter, §ion));
237cb8a476cSJames Wright PetscCall(PetscSectionSetFieldName(section, diff_filter_->field_prim_state, "Filtered Primitive State Variables"));
23888b07121SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_PRESSURE, "FilteredPressure"));
23988b07121SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_X, "FilteredVelocityX"));
24088b07121SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_Y, "FilteredVelocityY"));
24188b07121SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_VELOCITY_Z, "FilteredVelocityZ"));
24288b07121SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, DIFF_FILTER_TEMPERATURE, "FilteredTemperature"));
243cb8a476cSJames Wright PetscCall(PetscSectionSetFieldName(section, diff_filter_->field_velo_prod, "Filtered Velocity Products"));
2441ed803a2SJames Wright PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XX, "FilteredVelocitySquaredXX"));
2451ed803a2SJames Wright PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_YY, "FilteredVelocitySquaredYY"));
2461ed803a2SJames Wright PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_ZZ, "FilteredVelocitySquaredZZ"));
2471ed803a2SJames Wright PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_YZ, "FilteredVelocitySquaredYZ"));
2481ed803a2SJames Wright PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XZ, "FilteredVelocitySquaredXZ"));
2491ed803a2SJames Wright PetscCall(PetscSectionSetComponentName(section, 1, DIFF_FILTER_VELOCITY_SQUARED_XY, "FilteredVelocitySquaredXY"));
25025c92e8fSJames Wright }
25188b07121SJames Wright }
25288b07121SJames Wright
25388b07121SJames Wright PetscCall(PetscNew(&diff_filter_ctx));
254f5dc303cSJames Wright *diff_filter_ctx = (struct DifferentialFilterContext_){
255f5dc303cSJames Wright .grid_based_width = false,
256f5dc303cSJames Wright .width_scaling = {1, 1, 1},
257f5dc303cSJames Wright .kernel_scaling = 0.1,
258f5dc303cSJames Wright .damping_function = DIFF_FILTER_DAMP_NONE,
259f5dc303cSJames Wright .friction_length = 0,
260f5dc303cSJames Wright .damping_constant = 25,
261f5dc303cSJames Wright };
26288b07121SJames Wright
26388b07121SJames Wright PetscOptionsBegin(comm, NULL, "Differential Filtering Options", NULL);
26488b07121SJames Wright PetscInt narray = 3;
26588b07121SJames Wright PetscCall(PetscOptionsBool("-diff_filter_grid_based_width", "Use filter width based on the grid size", NULL, diff_filter_ctx->grid_based_width,
26688b07121SJames Wright (PetscBool *)&diff_filter_ctx->grid_based_width, NULL));
26788b07121SJames Wright PetscCall(PetscOptionsRealArray("-diff_filter_width_scaling", "Anisotropic scaling of filter width tensor", NULL, diff_filter_ctx->width_scaling,
26888b07121SJames Wright &narray, NULL));
26988b07121SJames Wright PetscCall(PetscOptionsReal("-diff_filter_kernel_scaling", "Scaling to make differential kernel size \"equivalent\" to other filter kernels", NULL,
27088b07121SJames Wright diff_filter_ctx->kernel_scaling, &diff_filter_ctx->kernel_scaling, NULL));
271682a9695SJames Wright PetscCall(PetscOptionsEnum("-diff_filter_wall_damping_function", "Damping function to use at the wall", NULL, DifferentialFilterDampingFunctions,
27214bd2a07SJames Wright (PetscEnum)diff_filter_ctx->damping_function, (PetscEnum *)&diff_filter_ctx->damping_function, NULL));
273682a9695SJames Wright PetscCall(PetscOptionsReal("-diff_filter_wall_damping_constant", "Contant for the wall-damping function", NULL, diff_filter_ctx->damping_constant,
274682a9695SJames Wright &diff_filter_ctx->damping_constant, NULL));
275682a9695SJames Wright PetscCall(PetscOptionsReal("-diff_filter_friction_length", "Friction length associated with the flow, \\delta_\\nu. For wall-damping functions",
276682a9695SJames Wright NULL, diff_filter_ctx->friction_length, &diff_filter_ctx->friction_length, NULL));
27788b07121SJames Wright PetscOptionsEnd();
27888b07121SJames Wright
2790c373b74SJames Wright Units units = honee->units;
28088b07121SJames Wright for (int i = 0; i < 3; i++) diff_filter_ctx->width_scaling[i] *= units->meter;
28188b07121SJames Wright diff_filter_ctx->kernel_scaling *= units->meter;
282682a9695SJames Wright diff_filter_ctx->friction_length *= units->meter;
28388b07121SJames Wright
28488b07121SJames Wright // -- Create QFContext
285cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx));
286cde3d787SJames Wright diff_filter_ctx->newt_ctx = *newt_ctx;
287cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx));
28888b07121SJames Wright
289b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(ceed, &diff_filter_qfctx));
290b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(diff_filter_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*diff_filter_ctx), diff_filter_ctx));
291b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(diff_filter_qfctx, CEED_MEM_HOST, FreeContextPetsc));
2929eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(diff_filter_qfctx, "filter width scaling",
2939eadbee4SJames Wright offsetof(struct DifferentialFilterContext_, width_scaling),
2949eadbee4SJames Wright sizeof(diff_filter_ctx->width_scaling) / sizeof(diff_filter_ctx->width_scaling[0]),
2959eadbee4SJames Wright "Filter width scaling"));
29688b07121SJames Wright
29788b07121SJames Wright // -- Setup Operators
298cb8a476cSJames Wright PetscCall(DifferentialFilterCreateOperators(honee, diff_filter_qfctx, diff_filter_));
299cb8a476cSJames Wright *diff_filter = diff_filter_;
30088b07121SJames Wright
301b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&diff_filter_qfctx));
302d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
30388b07121SJames Wright }
30488b07121SJames Wright
30588b07121SJames Wright // @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)306cb8a476cSJames Wright PetscErrorCode DifferentialFilterApply(Honee honee, DiffFilterData diff_filter, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution) {
30799dcded5SJames Wright PetscObjectState X_loc_state;
30899dcded5SJames Wright Vec RHS;
30988b07121SJames Wright
31088b07121SJames Wright PetscFunctionBeginUser;
311ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_DifferentialFilter, Q, Filtered_Solution, 0, 0));
31299dcded5SJames Wright PetscCall(DMGetNamedGlobalVector(diff_filter->dm_filter, "RHS", &RHS));
3130c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, diff_filter->op_rhs_ctx->X_loc, solution_time));
314ed6ba5acSJeremy L Thompson PetscCall(VecGetState(diff_filter->op_rhs_ctx->X_loc, &X_loc_state));
31599dcded5SJames Wright if (X_loc_state != diff_filter->X_loc_state) {
31699dcded5SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, RHS, diff_filter->op_rhs_ctx));
317ed6ba5acSJeremy L Thompson PetscCall(VecGetState(diff_filter->op_rhs_ctx->X_loc, &X_loc_state));
31899dcded5SJames Wright diff_filter->X_loc_state = X_loc_state;
31999dcded5SJames Wright }
32099dcded5SJames Wright PetscCall(VecViewFromOptions(RHS, NULL, "-diff_filter_rhs_view"));
32188b07121SJames Wright
32299dcded5SJames Wright PetscCall(KSPSolve(diff_filter->ksp, RHS, Filtered_Solution));
32399dcded5SJames Wright PetscCall(DMRestoreNamedGlobalVector(diff_filter->dm_filter, "RHS", &RHS));
324ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_DifferentialFilter, Q, Filtered_Solution, 0, 0));
325d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
32688b07121SJames Wright }
32788b07121SJames Wright
TSMonitor_DifferentialFilterSetup(TS ts,PetscViewerAndFormat * ctx)328cb8a476cSJames Wright PetscErrorCode TSMonitor_DifferentialFilterSetup(TS ts, PetscViewerAndFormat *ctx) {
329cb8a476cSJames Wright Honee honee;
330cb8a476cSJames Wright DiffFilterData diff_filter;
331cb8a476cSJames Wright
332cb8a476cSJames Wright PetscFunctionBeginUser;
333cb8a476cSJames Wright PetscCall(TSGetApplicationContext(ts, &honee));
334cb8a476cSJames Wright PetscCall(DifferentialFilterSetup(honee, &diff_filter));
335cb8a476cSJames Wright ctx->data = diff_filter;
336cb8a476cSJames Wright ctx->data_destroy = (PetscCtxDestroyFn *)DifferentialFilterDataDestroy;
337cb8a476cSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
338cb8a476cSJames Wright }
339cb8a476cSJames Wright
34088b07121SJames Wright // @brief TSMonitor for just applying differential filtering to the simulation
34188b07121SJames Wright // This runs every time step and is primarily for testing purposes
TSMonitor_DifferentialFilter(TS ts,PetscInt steps,PetscReal solution_time,Vec Q,PetscViewerAndFormat * ctx)342cb8a476cSJames Wright PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
343cb8a476cSJames Wright Honee honee;
344cb8a476cSJames Wright DiffFilterData diff_filter = (DiffFilterData)ctx->data;
34588b07121SJames Wright Vec Filtered_Field;
34688b07121SJames Wright
34788b07121SJames Wright PetscFunctionBeginUser;
348cb8a476cSJames Wright PetscCall(TSGetApplicationContext(ts, &honee));
34988b07121SJames Wright PetscCall(DMGetGlobalVector(diff_filter->dm_filter, &Filtered_Field));
35088b07121SJames Wright
351cb8a476cSJames Wright PetscCall(DifferentialFilterApply(honee, diff_filter, solution_time, Q, Filtered_Field));
35288b07121SJames Wright PetscCall(VecViewFromOptions(Filtered_Field, NULL, "-diff_filter_view"));
3530c373b74SJames Wright if (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER) PetscCall(RegressionTest(honee->app_ctx, Filtered_Field));
3541123f253SJames Wright
35588b07121SJames Wright PetscCall(DMRestoreGlobalVector(diff_filter->dm_filter, &Filtered_Field));
356d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
35788b07121SJames Wright }
35888b07121SJames Wright
DifferentialFilterDataDestroy(DiffFilterData * diff_filter)359cb8a476cSJames Wright PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData *diff_filter) {
360cb8a476cSJames Wright DiffFilterData diff_filter_ = *diff_filter;
361cb8a476cSJames Wright
36288b07121SJames Wright PetscFunctionBeginUser;
363cb8a476cSJames Wright if (!diff_filter_) PetscFunctionReturn(PETSC_SUCCESS);
364cb8a476cSJames Wright PetscCall(OperatorApplyContextDestroy(diff_filter_->op_rhs_ctx));
365cb8a476cSJames Wright PetscCall(DMDestroy(&diff_filter_->dm_filter));
366cb8a476cSJames Wright PetscCall(KSPDestroy(&diff_filter_->ksp));
367cb8a476cSJames Wright PetscCall(PetscFree(diff_filter_->num_field_components));
368cb8a476cSJames Wright PetscCall(PetscFree(diff_filter_));
369cb8a476cSJames Wright *diff_filter = NULL;
370d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
37188b07121SJames Wright }
37225c92e8fSJames Wright
DifferentialFilterMmsICSetup(Honee honee)3738d78d7c8SJames Wright PetscErrorCode DifferentialFilterMmsICSetup(Honee honee) {
37425c92e8fSJames Wright PetscFunctionBeginUser;
375f5dc303cSJames Wright honee->problem_data->ics = (HoneeQFSpec){
376f5dc303cSJames Wright .qf_func_ptr = DifferentialFilter_MMS_IC,
377f5dc303cSJames Wright .qf_loc = DifferentialFilter_MMS_IC_loc,
378f5dc303cSJames Wright .qfctx = honee->problem_data->ics.qfctx,
379f5dc303cSJames Wright };
380d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
38125c92e8fSJames Wright }
382