Line data Source code
1 0 : // Distributed under the MIT License.
2 : // See LICENSE.txt for details.
3 :
4 : #pragma once
5 :
6 : #include <array>
7 : #include <cstddef>
8 : #include <optional>
9 : #include <type_traits>
10 :
11 : #include "DataStructures/DataBox/DataBox.hpp"
12 : #include "DataStructures/DataBox/PrefixHelpers.hpp"
13 : #include "DataStructures/DataBox/Prefixes.hpp"
14 : #include "DataStructures/DataVector.hpp"
15 : #include "DataStructures/Tensor/Tensor.hpp"
16 : #include "DataStructures/Variables.hpp"
17 : #include "Domain/CoordinateMaps/Tags.hpp"
18 : #include "Domain/Structure/Element.hpp"
19 : #include "Domain/Tags.hpp"
20 : #include "Evolution/BoundaryCorrectionTags.hpp"
21 : #include "Evolution/DgSubcell/CartesianFluxDivergence.hpp"
22 : #include "Evolution/DgSubcell/ComputeBoundaryTerms.hpp"
23 : #include "Evolution/DgSubcell/CorrectPackagedData.hpp"
24 : #include "Evolution/DgSubcell/Tags/Coordinates.hpp"
25 : #include "Evolution/DgSubcell/Tags/Jacobians.hpp"
26 : #include "Evolution/DgSubcell/Tags/Mesh.hpp"
27 : #include "Evolution/DgSubcell/Tags/OnSubcellFaces.hpp"
28 : #include "Evolution/DiscontinuousGalerkin/Actions/PackageDataImpl.hpp"
29 : #include "Evolution/DiscontinuousGalerkin/MortarTags.hpp"
30 : #include "Evolution/Systems/ScalarAdvection/FiniteDifference/Reconstructor.hpp"
31 : #include "Evolution/Systems/ScalarAdvection/FiniteDifference/Tags.hpp"
32 : #include "Evolution/Systems/ScalarAdvection/Fluxes.hpp"
33 : #include "Evolution/Systems/ScalarAdvection/Subcell/ComputeFluxes.hpp"
34 : #include "Evolution/Systems/ScalarAdvection/System.hpp"
35 : #include "NumericalAlgorithms/Spectral/Mesh.hpp"
36 : #include "Utilities/CallWithDynamicType.hpp"
37 : #include "Utilities/ErrorHandling/Assert.hpp"
38 : #include "Utilities/Gsl.hpp"
39 : #include "Utilities/MakeArray.hpp"
40 :
41 : namespace ScalarAdvection::subcell {
42 : /*!
43 : * \brief Compute the time derivative on the subcell grid using FD
44 : * reconstruction.
45 : *
46 : * The code makes the following unchecked assumptions:
47 : * - Assumes Cartesian coordinates with a diagonal Jacobian matrix
48 : * from the logical to the inertial frame
49 : * - Assumes the mesh is not moving (grid and inertial frame are the same)
50 : */
51 : template <size_t Dim>
52 1 : struct TimeDerivative {
53 : template <typename DbTagsList>
54 0 : static void apply(const gsl::not_null<db::DataBox<DbTagsList>*> box) {
55 : ASSERT((db::get<::domain::CoordinateMaps::Tags::CoordinateMap<
56 : Dim, Frame::Grid, Frame::Inertial>>(*box))
57 : .is_identity(),
58 : "Do not yet support moving mesh with DG-subcell.");
59 : // subcell is currently not supported for external boundary elements
60 : const Element<Dim>& element = db::get<domain::Tags::Element<Dim>>(*box);
61 : ASSERT(element.external_boundaries().size() == 0,
62 : "Can't have external boundaries right now with subcell. ElementID "
63 : << element.id());
64 :
65 : using evolved_vars_tags = typename System<Dim>::variables_tag::tags_list;
66 : using fluxes_tags = typename Fluxes<Dim>::return_tags;
67 :
68 : // The copy of Mesh is intentional to avoid a GCC-7 internal compiler error.
69 : const Mesh<Dim> subcell_mesh =
70 : db::get<evolution::dg::subcell::Tags::Mesh<Dim>>(*box);
71 : ASSERT(
72 : subcell_mesh == Mesh<Dim>(subcell_mesh.extents(0),
73 : subcell_mesh.basis(0),
74 : subcell_mesh.quadrature(0)),
75 : "The subcell/FD mesh must be isotropic for the FD time derivative but "
76 : "got "
77 : << subcell_mesh);
78 :
79 : const size_t num_reconstructed_pts =
80 : (subcell_mesh.extents(0) + 1) *
81 : subcell_mesh.extents().slice_away(0).product();
82 :
83 : const ScalarAdvection::fd::Reconstructor<Dim>& recons =
84 : db::get<ScalarAdvection::fd::Tags::Reconstructor<Dim>>(*box);
85 :
86 : const auto& boundary_correction =
87 : db::get<evolution::Tags::BoundaryCorrection<System<Dim>>>(*box);
88 : using derived_boundary_corrections =
89 : typename std::decay_t<decltype(boundary_correction)>::creatable_classes;
90 :
91 : // Variables to store the boundary correction terms on FD subinterfaces
92 : std::array<Variables<evolved_vars_tags>, Dim> fd_boundary_corrections{};
93 :
94 : // package the data and compute the boundary correction
95 : tmpl::for_each<derived_boundary_corrections>([&boundary_correction,
96 : &fd_boundary_corrections,
97 : &box, &element,
98 : &num_reconstructed_pts,
99 : &recons, &subcell_mesh](
100 : auto
101 : derived_correction_v) {
102 : using derived_correction =
103 : tmpl::type_from<decltype(derived_correction_v)>;
104 : if (typeid(boundary_correction) == typeid(derived_correction)) {
105 : using dg_package_field_tags =
106 : typename derived_correction::dg_package_field_tags;
107 : using dg_package_data_temporary_tags =
108 : typename derived_correction::dg_package_data_temporary_tags;
109 : using dg_package_data_argument_tags =
110 : tmpl::append<evolved_vars_tags, fluxes_tags,
111 : dg_package_data_temporary_tags>;
112 :
113 : // Variables that need to be reconstructed for dg_package_data()
114 : auto package_data_argvars_lower_face = make_array<Dim>(
115 : Variables<dg_package_data_argument_tags>(num_reconstructed_pts));
116 : auto package_data_argvars_upper_face = make_array<Dim>(
117 : Variables<dg_package_data_argument_tags>(num_reconstructed_pts));
118 :
119 : // Reconstruct the fields on interfaces
120 : call_with_dynamic_type<void, typename ScalarAdvection::fd::
121 : Reconstructor<Dim>::creatable_classes>(
122 : &recons,
123 : [&box, &package_data_argvars_lower_face,
124 : &package_data_argvars_upper_face](const auto& reconstructor) {
125 : db::apply<typename std::decay_t<decltype(
126 : *reconstructor)>::reconstruction_argument_tags>(
127 : [&package_data_argvars_lower_face,
128 : &package_data_argvars_upper_face,
129 : &reconstructor](const auto&... args) {
130 : reconstructor->reconstruct(
131 : make_not_null(&package_data_argvars_lower_face),
132 : make_not_null(&package_data_argvars_upper_face),
133 : args...);
134 : },
135 : *box);
136 : });
137 :
138 : // copy over the face values of the velocity field
139 : {
140 : using tag = ::ScalarAdvection::Tags::VelocityField<Dim>;
141 : const auto& velocity_on_face =
142 : db::get<evolution::dg::subcell::Tags::OnSubcellFaces<tag, Dim>>(
143 : *box);
144 : for (size_t d = 0; d < Dim; ++d) {
145 : get<tag>(gsl::at(package_data_argvars_lower_face, d)) =
146 : gsl::at(velocity_on_face, d);
147 : get<tag>(gsl::at(package_data_argvars_upper_face, d)) =
148 : gsl::at(velocity_on_face, d);
149 : }
150 : }
151 :
152 : // Variables to store packaged data. Allocate outside of loop to
153 : // reduce allocations
154 : Variables<dg_package_field_tags> upper_packaged_data{
155 : num_reconstructed_pts};
156 : Variables<dg_package_field_tags> lower_packaged_data{
157 : num_reconstructed_pts};
158 :
159 : for (size_t dim = 0; dim < Dim; ++dim) {
160 : auto& vars_upper_face = gsl::at(package_data_argvars_upper_face, dim);
161 : auto& vars_lower_face = gsl::at(package_data_argvars_lower_face, dim);
162 :
163 : // Compute fluxes on each faces
164 : ScalarAdvection::subcell::compute_fluxes<Dim>(
165 : make_not_null(&vars_upper_face));
166 : ScalarAdvection::subcell::compute_fluxes<Dim>(
167 : make_not_null(&vars_lower_face));
168 :
169 : // Note that we use the sign convention on the normal vectors to
170 : // be compatible with DG.
171 : tnsr::i<DataVector, Dim, Frame::Inertial> lower_outward_conormal{
172 : num_reconstructed_pts, 0.0};
173 : lower_outward_conormal.get(dim) = 1.0;
174 : tnsr::i<DataVector, Dim, Frame::Inertial> upper_outward_conormal{
175 : num_reconstructed_pts, 0.0};
176 : upper_outward_conormal.get(dim) = -1.0;
177 :
178 : // Compute the packaged data
179 : evolution::dg::Actions::detail::dg_package_data<System<Dim>>(
180 : make_not_null(&upper_packaged_data),
181 : dynamic_cast<const derived_correction&>(boundary_correction),
182 : vars_upper_face, upper_outward_conormal, {std::nullopt}, *box,
183 : typename derived_correction::dg_package_data_volume_tags{},
184 : dg_package_data_argument_tags{});
185 : evolution::dg::Actions::detail::dg_package_data<System<Dim>>(
186 : make_not_null(&lower_packaged_data),
187 : dynamic_cast<const derived_correction&>(boundary_correction),
188 : vars_lower_face, lower_outward_conormal, {std::nullopt}, *box,
189 : typename derived_correction::dg_package_data_volume_tags{},
190 : dg_package_data_argument_tags{});
191 :
192 : // Now need to check if any of our neighbors are doing DG, because
193 : // if so then we need to use whatever boundary data they sent
194 : // instead of what we computed locally.
195 : //
196 : // Note: We could check this beforehand to avoid the extra work of
197 : // reconstruction and flux computations at the boundaries.
198 : evolution::dg::subcell::correct_package_data<true>(
199 : make_not_null(&lower_packaged_data),
200 : make_not_null(&upper_packaged_data), dim, element, subcell_mesh,
201 : db::get<evolution::dg::Tags::MortarData<Dim>>(*box), 0);
202 :
203 : // Compute the corrections on the faces. We only need to compute this
204 : // once because we can just flip the normal vectors then
205 : gsl::at(fd_boundary_corrections, dim)
206 : .initialize(num_reconstructed_pts);
207 : evolution::dg::subcell::compute_boundary_terms(
208 : make_not_null(&gsl::at(fd_boundary_corrections, dim)),
209 : dynamic_cast<const derived_correction&>(boundary_correction),
210 : upper_packaged_data, lower_packaged_data);
211 : }
212 : }
213 : });
214 :
215 : std::array<double, Dim> one_over_delta_xi{};
216 : {
217 : const tnsr::I<DataVector, Dim, Frame::ElementLogical>&
218 : cell_centered_logical_coords =
219 : db::get<evolution::dg::subcell::Tags::Coordinates<
220 : Dim, Frame::ElementLogical>>(*box);
221 :
222 : for (size_t i = 0; i < Dim; ++i) {
223 : // Note: assumes isotropic extents
224 : gsl::at(one_over_delta_xi, i) =
225 : 1.0 / (get<0>(cell_centered_logical_coords)[1] -
226 : get<0>(cell_centered_logical_coords)[0]);
227 : }
228 : }
229 :
230 : // Now compute the actual time derivative
231 : using dt_variables_tag =
232 : db::add_tag_prefix<::Tags::dt, typename System<Dim>::variables_tag>;
233 : const size_t num_pts = subcell_mesh.number_of_grid_points();
234 : db::mutate<dt_variables_tag>(
235 : [&num_pts, &fd_boundary_corrections, &subcell_mesh, &one_over_delta_xi](
236 : const auto dt_vars_ptr,
237 : const InverseJacobian<DataVector, Dim, Frame::ElementLogical,
238 : Frame::Grid>&
239 : cell_centered_logical_to_grid_inv_jacobian) {
240 : dt_vars_ptr->initialize(num_pts, 0.0);
241 : auto& dt_u =
242 : get<::Tags::dt<::ScalarAdvection::Tags::U>>(*dt_vars_ptr);
243 :
244 : for (size_t dim = 0; dim < Dim; ++dim) {
245 : Scalar<DataVector>& u_correction = get<::ScalarAdvection ::Tags::U>(
246 : gsl::at(fd_boundary_corrections, dim));
247 : evolution::dg::subcell::add_cartesian_flux_divergence(
248 : make_not_null(&get(dt_u)), gsl::at(one_over_delta_xi, dim),
249 : cell_centered_logical_to_grid_inv_jacobian.get(dim, dim),
250 : get(u_correction), subcell_mesh.extents(), dim);
251 : }
252 : },
253 : box,
254 : db::get<evolution::dg::subcell::fd::Tags::InverseJacobianLogicalToGrid<
255 : Dim>>(*box));
256 : }
257 : };
258 : } // namespace ScalarAdvection::subcell
|