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