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