NimbleSM
NimbleSM is a solid mechanics simulation code for dynamic systems
Loading...
Searching...
No Matches
nimble.mpi.reduction.h
Go to the documentation of this file.
1/*
2//@HEADER
3// ************************************************************************
4//
5// NimbleSM
6// Copyright 2018
7// National Technology & Engineering Solutions of Sandia, LLC (NTESS)
8//
9// Under the terms of Contract DE-NA0003525 with NTESS, the U.S. Government
10// retains certain rights in this software.
11//
12// Redistribution and use in source and binary forms, with or without
13// modification, are permitted provided that the following conditions are
14// met:
15//
16// 1. Redistributions of source code must retain the above copyright
17// notice, this list of conditions and the following disclaimer.
18//
19// 2. Redistributions in binary form must reproduce the above copyright
20// notice, this list of conditions and the following disclaimer in the
21// documentation and/or other materials provided with the distribution.
22//
23// 3. Neither the name of the Corporation nor the names of the
24// contributors may be used to endorse or promote products derived from
25// this software without specific prior written permission.
26//
27// THIS SOFTWARE IS PROVIDED BY NTESS "AS IS" AND ANY EXPRESS OR IMPLIED
28// WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF
29// MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN
30// NO EVENT SHALL NTESS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
31// SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
32// TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact David Littlewood (djlittl@sandia.gov)
39//
40// ************************************************************************
41//@HEADER
42*/
43
44#ifndef NIMBLE_MPI_REDUCTION_H
45#define NIMBLE_MPI_REDUCTION_H
46
47#ifdef NIMBLE_HAVE_MPI
48#include <mpi.h>
49
50#include <algorithm>
51#include <exception>
52#include <iostream>
53#include <memory>
54#include <numeric>
55#include <random>
56#include <set>
57#include <stdexcept>
58#include <type_traits>
59#include <unordered_map>
60#include <unordered_set>
61#include <vector>
62
65#include "nimble.quanta.h"
66
67namespace nimble {
68namespace reduction {
69struct ReductionInfo
70{
71 template <class... Args>
72 using vector = std::vector<Args...>;
73 template <class... Args>
74 using hashmap = std::unordered_map<Args...>;
75 std::vector<ReductionClique_t> cliques;
76 std::vector<int> unfinished;
77
78 ReductionInfo(
79 const vector<int>& clique_colors,
80 const vector<int>& clique_ids,
81 const vector<int>& clique_assignment,
82 const vector<int>& global_ids,
83 const mpicontext& context)
84 {
85 hashmap<int, vector<int>> clique_index_assignment;
86 clique_index_assignment.reserve(clique_ids.size());
87 int num_ranks = context.get_size();
88
89 hashmap<int, int> counts;
90 counts.reserve(clique_ids.size());
91
92 for (int clique_id : clique_assignment)
93 if (clique_id > num_ranks) counts[clique_id] += 1;
94
95 for (int clique_id : clique_ids) clique_index_assignment[clique_id].reserve(counts[clique_id]);
96
97 for (int index = 0, max = clique_assignment.size(); index < max; ++index) {
98 int clique_id = clique_assignment[index];
99 if (clique_id > num_ranks) { clique_index_assignment[clique_id].emplace_back(index); }
100 }
101
102 std::vector<MPI_Comm> comms{};
103 comms.reserve(clique_ids.size());
104
105 {
106 for (int color : clique_colors) {
107 MPI_Comm comm = context.split_by_color(color);
108 if (comm != MPI_COMM_NULL) { comms.push_back(comm); }
109 }
110 }
111
112 if (comms.size() != clique_ids.size()) NIMBLE_ABORT("**** Error, comms.size() != clique_ids.size().");
113
114 for (size_t i = 0; i < clique_ids.size(); ++i) {
115 int clique_id = clique_ids[i];
116 MPI_Comm comm = comms[i];
117 auto& index_list = clique_index_assignment[clique_id];
118 std::sort(index_list.begin(), index_list.end(), [&](int a, int b) { return global_ids[a] < global_ids[b]; });
119 this->cliques.emplace_back(std::move(index_list), index_list.size() * 3, comm);
120 }
121 }
122
123 ReductionInfo() {}
124 ReductionInfo(ReductionInfo&& ri) = default;
125 ReductionInfo&
126 operator=(ReductionInfo&& ri) = default;
127 template <int field_size, class Lookup>
128 void
129 Reduce(Lookup&& data)
130 {
131 for (auto& clique : cliques) clique.asyncreduce_initialize<field_size>(data);
132
133 unfinished.resize(cliques.size());
134 std::iota(unfinished.begin(), unfinished.end(), 0);
135
136 while (!unfinished.empty()) {
137 int increment = 0;
138 for (size_t i = 0; i < unfinished.size(); i += increment) {
139 bool reduceFinished = cliques[unfinished[i]].asyncreduce_finalize<field_size>(data);
140
141 increment = !reduceFinished;
142
143 if (reduceFinished) {
144 unfinished[i] = unfinished.back();
145 unfinished.pop_back();
146 }
147 }
148 }
149 }
150 void
151 GetAllIndices(std::vector<int>& indices, std::vector<int>& min_rank_containing_index)
152 {
153 for (auto& clique : cliques) {
154 int min_mpi_comm_world_rank = clique.GetMPICommWorldRanks()[0];
155 int num_indices = clique.GetNumIndices();
156 int const* clique_indices = clique.GetIndices();
157 for (int i = 0; i < num_indices; i++) {
158 indices.push_back(clique_indices[i]);
159 min_rank_containing_index.push_back(min_mpi_comm_world_rank);
160 }
161 }
162 }
163 void
164 PerformReduction(double* data, int field_size)
165 {
166 switch (field_size) {
167 case 1: Reduce<1>(data); break;
168 case 2: Reduce<2>(data); break;
169 case 3: Reduce<3>(data); break;
170 default:
171 std::string fs = std::to_string(field_size);
172 // std::cout << "Bad field size" << std::endl;
173 throw std::invalid_argument("Bad field size of " + fs);
174 }
175 }
176 template <class Lookup>
177 void
178 PerformReduction(Lookup& lookup, int field_size)
179 {
180 switch (field_size) {
181 case 1: Reduce<1>(lookup); break;
182 case 2: Reduce<2>(lookup); break;
183 case 3: Reduce<3>(lookup); break;
184 default:
185 std::string fs = std::to_string(field_size);
186 // std::cout << "Bad field size" << std::endl;
187 throw std::invalid_argument("Bad field size of " + fs);
188 }
189 }
190};
191
192template <class list_of_lists_t, class F>
193auto
194fill_clique_lookup(list_of_lists_t& ids_by_rank, F&& clique_lookup) ->
195 typename std::decay<quanta::transformed_iterated_t<F, quanta::elem_t<list_of_lists_t>>>::type
196{
197 typedef decltype(clique_lookup) lookup_t;
198 typedef quanta::elem_t<list_of_lists_t> inner_list_t;
199 typedef typename std::decay<quanta::transformed_iterated_t<lookup_t, inner_list_t>>::type clique_t;
200
201 std::unordered_map<clique_t, clique_t> remapped_cliques{};
202
203 const clique_t zero{};
204 clique_t rank_plus_one{};
205
206 auto generate_clique_id = quanta::make_counter<clique_t>(quanta::len(ids_by_rank) + 1);
207
208 for (auto& id_list : ids_by_rank) {
209 remapped_cliques.clear();
210 rank_plus_one++;
211
212 for (auto& id : id_list) {
213 auto& current_clique = clique_lookup(id);
214
215 if (current_clique == zero) {
216 current_clique = rank_plus_one;
217 } else {
218 auto remapped_clique_iter = remapped_cliques.find(current_clique);
219 if (remapped_clique_iter != remapped_cliques.end()) {
220 current_clique = remapped_clique_iter->second;
221 } else {
222 auto new_clique = generate_clique_id();
223 remapped_cliques[current_clique] = new_clique;
224 current_clique = new_clique;
225 }
226 }
227 }
228 }
229 return generate_clique_id.get_count();
230}
231
232ReductionInfo*
233GenerateReductionInfo(const std::vector<int>& raw_global_ids, const mpicontext& context);
234
235} // namespace reduction
236} // namespace nimble
237#endif
238
239#endif // NIMBLE_MPI_REDUCTION_H
Definition kokkos_contact_manager.h:49
#define NIMBLE_ABORT(...)
Definition nimble_macros.h:87