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