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