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