SpECTRE Documentation Coverage Report
Current view: top level - ParallelAlgorithms/LinearSolver/Multigrid/Actions - RestrictFields.hpp Hit Total Coverage
Commit: de0084593a37be7727c6c4da0be2184b0f8d9ed4 Lines: 4 8 50.0 %
Date: 2025-11-04 23:26:01
Legend: Lines: hit not hit

          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 <map>
       8             : #include <optional>
       9             : #include <tuple>
      10             : #include <type_traits>
      11             : #include <utility>
      12             : 
      13             : #include "DataStructures/ApplyMatrices.hpp"
      14             : #include "DataStructures/DataBox/DataBox.hpp"
      15             : #include "DataStructures/DataBox/PrefixHelpers.hpp"
      16             : #include "DataStructures/FixedHashMap.hpp"
      17             : #include "Domain/Structure/ChildSize.hpp"
      18             : #include "Domain/Structure/ElementId.hpp"
      19             : #include "Domain/Tags.hpp"
      20             : #include "IO/Logging/Tags.hpp"
      21             : #include "IO/Logging/Verbosity.hpp"
      22             : #include "IO/Observer/Tags.hpp"
      23             : #include "NumericalAlgorithms/Convergence/Tags.hpp"
      24             : #include "NumericalAlgorithms/Spectral/Projection.hpp"
      25             : #include "Parallel/AlgorithmExecution.hpp"
      26             : #include "Parallel/GlobalCache.hpp"
      27             : #include "Parallel/InboxInserters.hpp"
      28             : #include "Parallel/Invoke.hpp"
      29             : #include "Parallel/Printf/Printf.hpp"
      30             : #include "ParallelAlgorithms/Amr/Tags.hpp"
      31             : #include "ParallelAlgorithms/LinearSolver/Multigrid/Tags.hpp"
      32             : #include "ParallelAlgorithms/LinearSolver/Tags.hpp"
      33             : #include "Utilities/ConstantExpressions.hpp"
      34             : #include "Utilities/ErrorHandling/Assert.hpp"
      35             : #include "Utilities/GetOutput.hpp"
      36             : #include "Utilities/Gsl.hpp"
      37             : #include "Utilities/PrettyType.hpp"
      38             : #include "Utilities/TMPL.hpp"
      39             : #include "Utilities/TaggedTuple.hpp"
      40             : 
      41             : /// \cond
      42             : template <size_t Dim>
      43             : struct ElementId;
      44             : /// \endcond
      45             : 
      46           1 : namespace LinearSolver::multigrid {
      47             : 
      48             : template <size_t Dim, typename ReceiveTags>
      49           0 : struct DataFromChildrenInboxTag
      50             :     : public Parallel::InboxInserters::Map<
      51             :           DataFromChildrenInboxTag<Dim, ReceiveTags>> {
      52           0 :   using temporal_id = size_t;
      53           0 :   using type =
      54             :       std::map<temporal_id,
      55             :                FixedHashMap<two_to_the(Dim), ElementId<Dim>,
      56             :                             tuples::tagged_tuple_from_typelist<ReceiveTags>,
      57             :                             boost::hash<ElementId<Dim>>>>;
      58             : };
      59             : 
      60             : /// Actions related to the Multigrid linear solver
      61           1 : namespace Actions {
      62             : /*!
      63             :  * \brief Communicate and project the `FieldsTags` to the next-coarser grid
      64             :  * in the multigrid hierarchy
      65             :  *
      66             :  * \tparam FieldsTags These tags will be communicated and projected. They can
      67             :  * hold any type that works with `::apply_matrices` and supports addition, e.g.
      68             :  * `Variables`.
      69             :  * \tparam OptionsGroup The option group identifying the multigrid solver
      70             :  * \tparam FieldsAreMassiveTag A boolean tag in the DataBox that indicates
      71             :  * whether or not the `FieldsTags` have already been multiplied by the mass
      72             :  * matrix. This setting influences the way the fields are projected. In
      73             :  * particular, the mass matrix already includes a Jacobian factor, so the
      74             :  * difference in size between the parent and the child element is already
      75             :  * accounted for.
      76             :  * \tparam ReceiveTags The projected fields will be stored in these tags
      77             :  * (default: `FieldsTags`).
      78             :  */
      79             : template <typename FieldsTags, typename OptionsGroup,
      80             :           typename FieldsAreMassiveTag, typename ReceiveTags = FieldsTags>
      81           1 : struct SendFieldsToCoarserGrid;
      82             : 
      83             : /// \cond
      84             : template <typename... FieldsTags, typename OptionsGroup,
      85             :           typename FieldsAreMassiveTag, typename... ReceiveTags>
      86             : struct SendFieldsToCoarserGrid<tmpl::list<FieldsTags...>, OptionsGroup,
      87             :                                FieldsAreMassiveTag,
      88             :                                tmpl::list<ReceiveTags...>> {
      89             :   using const_global_cache_tags =
      90             :       tmpl::list<logging::Tags::Verbosity<OptionsGroup>>;
      91             : 
      92             :   template <typename DbTagsList, typename... InboxTags, typename Metavariables,
      93             :             size_t Dim, typename ActionList, typename ParallelComponent>
      94             :   static Parallel::iterable_action_return_t apply(
      95             :       db::DataBox<DbTagsList>& box,
      96             :       const tuples::TaggedTuple<InboxTags...>& /*inboxes*/,
      97             :       Parallel::GlobalCache<Metavariables>& cache,
      98             :       const ElementId<Dim>& element_id, const ActionList /*meta*/,
      99             :       const ParallelComponent* const /*meta*/) {
     100             :     // Skip restriction on coarsest level
     101             :     const auto& parent_id = db::get<amr::Tags::ParentId<Dim>>(box);
     102             :     if (not parent_id.has_value()) {
     103             :       return {Parallel::AlgorithmExecution::Continue, std::nullopt};
     104             :     }
     105             : 
     106             :     const size_t iteration_id =
     107             :         db::get<Convergence::Tags::IterationId<OptionsGroup>>(box);
     108             :     if (UNLIKELY(db::get<logging::Tags::Verbosity<OptionsGroup>>(box) >=
     109             :                  ::Verbosity::Debug)) {
     110             :       Parallel::printf("%s %s(%zu): Send fields to coarser grid %s\n",
     111             :                        element_id, pretty_type::name<OptionsGroup>(),
     112             :                        iteration_id, parent_id);
     113             :     }
     114             : 
     115             :     // Restrict the fields to the coarser (parent) grid.
     116             :     // We restrict before sending the data so the restriction operation is
     117             :     // parellelized. The parent only needs to sum up all child contributions.
     118             :     const auto& mesh = db::get<domain::Tags::Mesh<Dim>>(box);
     119             :     const auto& parent_mesh = db::get<amr::Tags::ParentMesh<Dim>>(box);
     120             :     ASSERT(
     121             :         parent_mesh.has_value(),
     122             :         "Should have a parent mesh, because a parent ID is set. This element: "
     123             :             << element_id << ", parent element: " << *parent_id);
     124             :     const auto child_size =
     125             :         domain::child_size(element_id.segment_ids(), parent_id->segment_ids());
     126             :     bool massive = false;
     127             :     if constexpr (not std::is_same_v<FieldsAreMassiveTag, void>) {
     128             :       massive = db::get<FieldsAreMassiveTag>(box);
     129             :     }
     130             :     tuples::TaggedTuple<ReceiveTags...> restricted_fields{};
     131             :     if (Spectral::needs_projection(mesh, *parent_mesh, child_size)) {
     132             :       const auto restriction_operator =
     133             :           Spectral::projection_matrix_child_to_parent(mesh, *parent_mesh,
     134             :                                                       child_size, massive);
     135             :       const auto restrict_fields = [&restricted_fields, &restriction_operator,
     136             :                                     &mesh](const auto receive_tag_v,
     137             :                                            const auto& fields) {
     138             :         using receive_tag = std::decay_t<decltype(receive_tag_v)>;
     139             :         get<receive_tag>(restricted_fields) = typename receive_tag::type(
     140             :             apply_matrices(restriction_operator, fields, mesh.extents()));
     141             :         return '0';
     142             :       };
     143             :       expand_pack(restrict_fields(ReceiveTags{}, db::get<FieldsTags>(box))...);
     144             :     } else {
     145             :       expand_pack(
     146             :           (get<ReceiveTags>(restricted_fields) =
     147             :                typename ReceiveTags::type(db::get<FieldsTags>(box)))...);
     148             :     }
     149             : 
     150             :     // Send restricted fields to the parent
     151             :     auto& receiver_proxy =
     152             :         Parallel::get_parallel_component<ParallelComponent>(cache);
     153             :     Parallel::receive_data<
     154             :         DataFromChildrenInboxTag<Dim, tmpl::list<ReceiveTags...>>>(
     155             :         receiver_proxy[*parent_id], iteration_id,
     156             :         std::make_pair(element_id, std::move(restricted_fields)));
     157             :     return {Parallel::AlgorithmExecution::Continue, std::nullopt};
     158             :   }
     159             : };
     160             : /// \endcond
     161             : 
     162             : /// Receive the `FieldsTags` communicated from the finer grid in the multigrid
     163             : /// hierarchy.
     164             : ///
     165             : /// \see LinearSolver::multigrid::Actions::SendFieldsToCoarserGrid
     166             : template <size_t Dim, typename FieldsTags, typename OptionsGroup,
     167             :           typename ReceiveTags = FieldsTags>
     168           1 : struct ReceiveFieldsFromFinerGrid;
     169             : 
     170             : /// \cond
     171             : template <size_t Dim, typename FieldsTags, typename OptionsGroup,
     172             :           typename... ReceiveTags>
     173             : struct ReceiveFieldsFromFinerGrid<Dim, FieldsTags, OptionsGroup,
     174             :                                   tmpl::list<ReceiveTags...>> {
     175             :   using inbox_tags =
     176             :       tmpl::list<DataFromChildrenInboxTag<Dim, tmpl::list<ReceiveTags...>>>;
     177             :   using const_global_cache_tags =
     178             :       tmpl::list<logging::Tags::Verbosity<OptionsGroup>>;
     179             : 
     180             :   template <typename DbTagsList, typename... InboxTags, typename Metavariables,
     181             :             typename ActionList, typename ParallelComponent>
     182             :   static Parallel::iterable_action_return_t apply(
     183             :       db::DataBox<DbTagsList>& box, tuples::TaggedTuple<InboxTags...>& inboxes,
     184             :       const Parallel::GlobalCache<Metavariables>& /*cache*/,
     185             :       const ElementId<Dim>& element_id, const ActionList /*meta*/,
     186             :       const ParallelComponent* const /*meta*/) {
     187             :     // Skip on finest grid
     188             :     const auto& child_ids = db::get<amr::Tags::ChildIds<Dim>>(box);
     189             :     const bool is_finest_grid = child_ids.empty();
     190             :     if (is_finest_grid) {
     191             :       return {Parallel::AlgorithmExecution::Continue, std::nullopt};
     192             :     }
     193             : 
     194             :     // Wait for data from finer grid
     195             :     const size_t iteration_id =
     196             :         db::get<Convergence::Tags::IterationId<OptionsGroup>>(box);
     197             :     auto& inbox =
     198             :         tuples::get<DataFromChildrenInboxTag<Dim, tmpl::list<ReceiveTags...>>>(
     199             :             inboxes);
     200             :     const auto received_this_iteration = inbox.find(iteration_id);
     201             :     if (received_this_iteration == inbox.end()) {
     202             :       if (UNLIKELY(db::get<logging::Tags::Verbosity<OptionsGroup>>(box) >=
     203             :                    ::Verbosity::Debug)) {
     204             :         Parallel::printf(
     205             :             "%s %s(%zu): Waiting for fine-grid data (still empty)\n",
     206             :             element_id, pretty_type::name<OptionsGroup>(), iteration_id);
     207             :       }
     208             :       return {Parallel::AlgorithmExecution::Retry, std::nullopt};
     209             :     }
     210             :     const auto& received_children_data = received_this_iteration->second;
     211             :     for (const auto& child_id : child_ids) {
     212             :       if (received_children_data.find(child_id) ==
     213             :           received_children_data.end()) {
     214             :         if (UNLIKELY(db::get<logging::Tags::Verbosity<OptionsGroup>>(box) >=
     215             :                      ::Verbosity::Debug)) {
     216             :           Parallel::printf(
     217             :               "%s %s(%zu): Waiting for fine-grid data from child %s\n",
     218             :               element_id, pretty_type::name<OptionsGroup>(), iteration_id,
     219             :               child_id);
     220             :         }
     221             :         return {Parallel::AlgorithmExecution::Retry, std::nullopt};
     222             :       }
     223             :     }
     224             :     auto children_data = std::move(inbox.extract(iteration_id).mapped());
     225             : 
     226             :     if (UNLIKELY(db::get<logging::Tags::Verbosity<OptionsGroup>>(box) >=
     227             :                  ::Verbosity::Debug)) {
     228             :       Parallel::printf("%s %s(%zu): Receive fields from finer grid\n",
     229             :                        element_id, pretty_type::name<OptionsGroup>(),
     230             :                        iteration_id);
     231             :     }
     232             : 
     233             :     // Assemble restricted data from children
     234             :     const auto assemble_children_data =
     235             :         [&children_data](const auto source, const auto receive_tag_v) {
     236             :           using receive_tag = std::decay_t<decltype(receive_tag_v)>;
     237             :           // Move the first child data directly into the buffer, then add the
     238             :           // data from the remaining children.
     239             :           auto child_id_and_data = children_data.begin();
     240             :           *source = std::move(get<receive_tag>(child_id_and_data->second));
     241             :           ++child_id_and_data;
     242             :           while (child_id_and_data != children_data.end()) {
     243             :             *source += get<receive_tag>(child_id_and_data->second);
     244             :             ++child_id_and_data;
     245             :           }
     246             :           return '0';
     247             :         };
     248             :     expand_pack(db::mutate<ReceiveTags>(assemble_children_data,
     249             :                                         make_not_null(&box), ReceiveTags{})...);
     250             : 
     251             :     return {Parallel::AlgorithmExecution::Continue, std::nullopt};
     252             :   }
     253             : };
     254             : /// \endcond
     255             : 
     256             : }  // namespace Actions
     257             : }  // namespace LinearSolver::multigrid

Generated by: LCOV version 1.14