Line data Source code
1 0 : // Distributed under the MIT License. 2 : // See LICENSE.txt for details. 3 : 4 : #pragma once 5 : 6 : #include <cstddef> 7 : #include <memory> 8 : #include <optional> 9 : #include <pup.h> 10 : #include <string> 11 : #include <unordered_map> 12 : #include <unordered_set> 13 : #include <vector> 14 : 15 : #include "DataStructures/Tensor/TypeAliases.hpp" 16 : #include "Utilities/Gsl.hpp" 17 : #include "Utilities/Serialization/CharmPupable.hpp" 18 : #include "Utilities/TMPL.hpp" 19 : 20 : /// \cond 21 : class DataVector; 22 : template <size_t Dim> 23 : class Mesh; 24 : template <typename TagsList> 25 : class Variables; 26 : /// \endcond 27 : 28 : namespace Filters { 29 : /*! 30 : * \ingroup DiscontinuousGalerkinGroup 31 : * \brief Abstract base class for spectral filters. 32 : * 33 : * Derived classes implement spectral filters that suppress high-frequency 34 : * modes in the modal expansion of a discontinuous Galerkin solution. Writing 35 : * the 1-D modal coefficients of a tensor component as \f$c_i\f$, a spectral 36 : * filter rescales them as 37 : * 38 : * \f{align*}{ 39 : * c_i \to \sigma\!\left(\frac{i}{N}\right) c_i, 40 : * \f} 41 : * 42 : * where \f$N\f$ is the basis degree (number of grid points per element per 43 : * dimension minus one) and \f$\sigma\f$ is a derived-class-specific damping 44 : * factor (for example, an exponential in `Filters::Hypercube`). For a 45 : * discussion of filtering see section 5.3 of \cite HesthavenWarburton. 46 : * 47 : * The interface is organized along two axes: 48 : * 49 : * - **Volume vs. boundary.** `apply_in_volume` filters the evolved 50 : * `Variables` on the element's `Dim`-dimensional mesh, while 51 : * `apply_on_boundary` filters data on a `Dim - 1` face mesh. The boundary 52 : * hook is intended for filtering the lifted boundary corrections produced 53 : * by the DG numerical flux before they enter the volume right-hand side. 54 : * 55 : * - **Substep vs. every-N-steps.** The driving action queries 56 : * `apply_volume_filter_on_substep()` (and the boundary analogue) to decide 57 : * whether to filter inside every Runge-Kutta substep, and 58 : * `apply_volume_filter_on_this_step(step_number)` (and the boundary 59 : * analogue) to gate filtering at coarser cadences such as every `N` full 60 : * steps. The two return values are independent so a derived class may pick 61 : * any combination. 62 : * 63 : * Derived classes additionally declare via `need_jacobians()` whether they 64 : * require the grid-to-inertial Jacobian and its inverse in `apply_in_volume` 65 : * / `apply_on_boundary`. When `need_jacobians()` returns `false`, the 66 : * driving action passes `std::nullopt` for both Jacobian arguments and skips 67 : * computing them. 68 : * 69 : * The set of domain blocks to which the filter applies is reported by 70 : * `blocks_to_filter()`; a `std::nullopt` return value means the filter 71 : * applies to every block in the domain. 72 : */ 73 : template <size_t Dim, typename TagList> 74 1 : class Filter : public PUP::able { 75 : public: 76 0 : Filter() = default; 77 0 : ~Filter() noexcept override = default; 78 : 79 0 : WRAPPED_PUPable_abstract(Filter); // NOLINT 80 0 : explicit Filter(CkMigrateMessage* m) : PUP::able(m) {} 81 : 82 : /// Return a heap-allocated deep copy of this filter. 83 1 : virtual std::unique_ptr<Filter<Dim, TagList>> get_clone() const = 0; 84 : 85 : /// Whether the volume filter should be applied inside every Runge-Kutta 86 : /// substep. 87 1 : virtual bool apply_volume_filter_on_substep() const = 0; 88 : 89 : /// Whether the volume filter should be applied at the given `step_number`. 90 1 : virtual bool apply_volume_filter_on_this_step(size_t step_number) const = 0; 91 : 92 : /// Whether the boundary filter should be applied inside every Runge-Kutta 93 : /// substep. 94 1 : virtual bool apply_boundary_filter_on_substep() const = 0; 95 : 96 : /// Whether the boundary filter should be applied at the given 97 : /// `step_number`. 98 1 : virtual bool apply_boundary_filter_on_this_step(size_t step_number) const = 0; 99 : 100 : /// \brief Whether the filter needs the grid-to-inertial Jacobian and its 101 : /// inverse. 102 : /// 103 : /// When `false`, the driving action passes `std::nullopt` for the Jacobian 104 : /// arguments of `apply_in_volume` and `apply_on_boundary` and avoids 105 : /// constructing them. 106 1 : virtual bool need_jacobians() const = 0; 107 : 108 : /// \brief Returns `true` if this filter can filter the `mesh`. 109 1 : virtual bool supports_mesh(const Mesh<Dim>& mesh) const = 0; 110 : 111 : /// \brief Returns `true` if `other` is the same concrete filter type and is 112 : /// equivalent to `*this`. 113 1 : virtual bool is_equal(const Filter& other) const = 0; 114 : 115 : /// \brief The set of block (or block group) names the filter applies to. 116 : /// 117 : /// `std::nullopt` means the filter applies to every block in the domain. 118 1 : virtual const std::optional<std::vector<size_t>>& blocks_to_filter() 119 : const = 0; 120 : 121 : /// \brief Used after construction to change blocks and groups labeled by 122 : /// strings to `size_t`s, which is what is used throughout the code. 123 1 : virtual void set_blocks_to_filter( 124 : const std::vector<std::string>& all_block_names, 125 : const std::unordered_map<std::string, std::unordered_set<std::string>>& 126 : block_groups) = 0; 127 : 128 : /*! 129 : * \brief Apply the filter in place to the volume `vars` on the 130 : * `Dim`-dimensional `mesh`. 131 : * 132 : * The grid-to-inertial Jacobian `jac_grid_to_inertial` and its inverse 133 : * `inv_jac_grid_to_inertial` are populated only when `need_jacobians()` 134 : * returns `true`; otherwise both arguments are `std::nullopt`. 135 : */ 136 1 : virtual void apply_in_volume( 137 : gsl::not_null<Variables<TagList>*> vars, const Mesh<Dim>& mesh, 138 : const std::optional< 139 : InverseJacobian<DataVector, Dim, Frame::Grid, Frame::Inertial>>& 140 : inv_jac_grid_to_inertial, 141 : const std::optional< 142 : Jacobian<DataVector, Dim, Frame::Grid, Frame::Inertial>>& 143 : jac_grid_to_inertial) const = 0; 144 : 145 : /*! 146 : * \brief Apply the filter in place to face `vars` on the `Dim - 1` face 147 : * `mesh`. 148 : * 149 : * Used to filter the lifted boundary corrections produced by the DG numerical 150 : * flux before they enter the volume right-hand side. The grid-to-inertial 151 : * Jacobian `jac_grid_to_inertial` and its inverse `inv_jac_grid_to_inertial` 152 : * are populated only when `need_jacobians()` returns `true`; otherwise both 153 : * arguments are `std::nullopt`. 154 : */ 155 1 : virtual void apply_on_boundary( 156 : gsl::not_null<Variables<TagList>*> vars, const Mesh<Dim - 1>& mesh, 157 : const std::optional< 158 : InverseJacobian<DataVector, Dim, Frame::Grid, Frame::Inertial>>& 159 : inv_jac_grid_to_inertial, 160 : const std::optional< 161 : Jacobian<DataVector, Dim, Frame::Grid, Frame::Inertial>>& 162 : jac_grid_to_inertial) const = 0; 163 : }; 164 : } // namespace Filters