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