ThunderEgg  1.0.0
Interface.h
Go to the documentation of this file.
1 /***************************************************************************
2  * ThunderEgg, a library for solvers on adaptively refined block-structured
3  * Cartesian grids.
4  *
5  * Copyright (c) 2020-2021 Scott Aiton
6  *
7  * This program is free software: you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation, either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program. If not, see <https://www.gnu.org/licenses/>.
19  ***************************************************************************/
20 
21 #ifndef THUNDEREGG_SCHUR_INTERFACE_H
22 #define THUNDEREGG_SCHUR_INTERFACE_H
23 
32 #include <bitset>
33 #include <map>
34 #include <mpi.h>
35 #include <set>
36 #include <vector>
37 namespace ThunderEgg {
38 namespace Schur {
46 template<int D>
47 class Interface : public Serializable
48 {
49 public:
53  int id = -1;
57  int local_index = -1;
61  int global_index = -1;
69  {
70  public:
82  std::shared_ptr<const PatchIfaceInfo<D>> piinfo;
83 
89  std::shared_ptr<PatchIfaceInfo<D>> getNonConstPiinfo()
90  {
91  return std::const_pointer_cast<PatchIfaceInfo<D>>(piinfo);
92  }
93 
102  : side(side)
103  , type(type)
104  , piinfo(piinfo)
105  {}
106  };
113  std::vector<SideTypePiinfo> patches;
114 
118  Interface() = default;
119 
125  explicit Interface(int id)
126  : id(id)
127  {}
128 
135  void insert(Side<D> s, std::shared_ptr<const PatchIfaceInfo<D>> piinfo)
136  {
137  NbrType nbr_type = piinfo->pinfo.getNbrType(s);
138  if (nbr_type == NbrType::Normal) {
139  patches.emplace_back(s, IfaceType<D>::Normal(), piinfo);
140  } else if (nbr_type == NbrType::Fine) {
141  auto info = piinfo->getFineIfaceInfo(s);
142  if (info->id == id) {
143  patches.emplace_back(s, IfaceType<D>::CoarseToCoarse(), piinfo);
144  } else {
146  if (info->fine_ids[o.getIndex()] == id) {
147  patches.emplace_back(s, IfaceType<D>::CoarseToFine(o), piinfo);
148  break;
149  }
150  }
151  }
152  } else if (nbr_type == NbrType::Coarse) {
153  auto info = piinfo->getCoarseIfaceInfo(s);
154  if (info->id == id) {
155  patches.emplace_back(s, IfaceType<D>::FineToFine(info->orth_on_coarse), piinfo);
156  } else {
157  patches.emplace_back(s, IfaceType<D>::FineToCoarse(info->orth_on_coarse), piinfo);
158  }
159  } else {
160  throw RuntimeError("Unsupported NbrType");
161  }
162  }
168  void merge(Interface<D> ifs)
169  {
170  for (auto patch : ifs.patches) {
171  patches.push_back(patch);
172  }
173  }
174  int serialize(char* buffer) const
175  {
176  BufferWriter writer(buffer);
177  writer << id;
178  writer << global_index;
179  int size = (int)patches.size();
180  writer << size;
181  for (auto patch : patches) {
182  writer << patch.side;
183  writer << patch.type;
184  writer << *patch.piinfo;
185  }
186  return writer.getPos();
187  }
188  int deserialize(char* buffer)
189  {
190  BufferReader reader(buffer);
191  reader >> id;
192  reader >> global_index;
193  int size = 0;
194  reader >> size;
195  for (int i = 0; i < size; i++) {
196  Side<D> s;
197  IfaceType<D> type;
198  auto piinfo = std::make_shared<PatchIfaceInfo<D>>();
199  reader >> s;
200  reader >> type;
201  reader >> *piinfo;
202  patches.emplace_back(s, type, piinfo);
203  }
204  return reader.getPos();
205  }
206 
207 private:
217  static void InsertPatchToInterface(
218  std::map<int, std::map<int, std::shared_ptr<Interface<D>>>>& rank_id_iface_map,
219  int rank,
220  int id,
221  Side<D> s,
222  std::shared_ptr<const PatchIfaceInfo<D>> piinfo)
223  {
224  std::shared_ptr<Interface<D>>& iface_ptr = rank_id_iface_map[rank][id];
225  if (iface_ptr == nullptr) {
226  iface_ptr.reset(new Interface<D>(id));
227  }
228  iface_ptr->insert(s, piinfo);
229  }
239  static void InsertInterfaceWithNormalNbr(
240  std::map<int, std::map<int, std::shared_ptr<Interface<D>>>>& rank_id_iface_map,
241  std::set<int>& incoming_procs,
242  std::shared_ptr<const PatchIfaceInfo<D>> piinfo,
243  Side<D> s)
244  {
245  int rank;
246  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
247  auto info = piinfo->getNormalIfaceInfo(s);
248  auto nbr_info = piinfo->pinfo.getNormalNbrInfo(s);
249  InsertPatchToInterface(rank_id_iface_map, info->rank, info->id, s, piinfo);
250  if (info->rank == piinfo->pinfo.rank && nbr_info.rank != piinfo->pinfo.rank) {
251  incoming_procs.insert(nbr_info.rank);
252  }
253  }
263  static void InsertInterfacesWithFineNbr(
264  std::map<int, std::map<int, std::shared_ptr<Interface<D>>>>& rank_id_iface_map,
265  std::set<int>& incoming_procs,
266  std::shared_ptr<const PatchIfaceInfo<D>> piinfo,
267  Side<D> s)
268  {
269  auto info = piinfo->getFineIfaceInfo(s);
270 
271  InsertPatchToInterface(rank_id_iface_map, info->rank, info->id, s, piinfo);
272 
273  for (size_t i = 0; i < Orthant<D - 1>::num_orthants; i++) {
274  InsertPatchToInterface(rank_id_iface_map, info->fine_ranks[i], info->fine_ids[i], s, piinfo);
275 
276  if (info->fine_ranks[i] != piinfo->pinfo.rank) {
277  incoming_procs.insert(info->fine_ranks[i]);
278  }
279  }
280  }
290  static void InsertInterfacesWithCoarseNbr(
291  std::map<int, std::map<int, std::shared_ptr<Interface<D>>>>& rank_id_iface_map,
292  std::set<int>& incoming_procs,
293  std::shared_ptr<const PatchIfaceInfo<D>> piinfo,
294  Side<D> s)
295  {
296  auto info = piinfo->getCoarseIfaceInfo(s);
297 
298  InsertPatchToInterface(rank_id_iface_map, info->rank, info->id, s, piinfo);
299  InsertPatchToInterface(rank_id_iface_map, info->coarse_rank, info->coarse_id, s, piinfo);
300 
301  if (info->coarse_rank != piinfo->pinfo.rank) {
302  incoming_procs.insert(info->coarse_rank);
303  }
304  }
305 
306 public:
318  std::vector<std::shared_ptr<const PatchIfaceInfo<D>>> piinfos,
319  std::map<int, std::map<int, std::shared_ptr<Interface<D>>>>& rank_id_iface_map,
320  std::vector<std::shared_ptr<PatchIfaceInfo<D>>>& off_proc_piinfos)
321  {
322  int rank;
323  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
324  rank_id_iface_map.clear();
325  std::set<int> incoming_procs;
326  for (auto piinfo : piinfos) {
327  for (Side<D> s : Side<D>::getValues()) {
328  if (piinfo->pinfo.hasNbr(s)) {
329  switch (piinfo->pinfo.getNbrType(s)) {
330  case NbrType::Normal:
331  InsertInterfaceWithNormalNbr(rank_id_iface_map, incoming_procs, piinfo, s);
332  break;
333  case NbrType::Fine:
334  InsertInterfacesWithFineNbr(rank_id_iface_map, incoming_procs, piinfo, s);
335  break;
336  case NbrType::Coarse:
337  InsertInterfacesWithCoarseNbr(rank_id_iface_map, incoming_procs, piinfo, s);
338  break;
339  default:
340  throw RuntimeError("Unsupported NbrType value");
341  }
342  }
343  }
344  }
345  // send info
346  std::deque<std::vector<char>> buffers;
347  std::vector<MPI_Request> send_requests;
348  for (auto& p : rank_id_iface_map) {
349  int dest = p.first;
350  if (dest != rank) {
351  int size = 0;
352  for (auto q : p.second) {
353  auto& iface = q.second;
354  size += iface->serialize(nullptr);
355  }
356  buffers.emplace_back(size);
357  BufferWriter writer(buffers.back().data());
358  for (auto q : p.second) {
359  auto& iface = q.second;
360  writer << *iface;
361  }
362  MPI_Request request;
363  MPI_Isend(buffers.back().data(), size, MPI_BYTE, dest, 0, MPI_COMM_WORLD, &request);
364  send_requests.push_back(request);
365  }
366  }
367  // recv info
368  std::map<int, std::shared_ptr<PatchIfaceInfo<D>>> id_to_off_proc_piinfo_map;
369 
370  MPI_Barrier(MPI_COMM_WORLD);
371 
372  size_t num_incoming = incoming_procs.size();
373  for (size_t i = 0; i < num_incoming; i++) {
374  MPI_Status status;
375  MPI_Probe(MPI_ANY_SOURCE, 0, MPI_COMM_WORLD, &status);
376  int size;
377  MPI_Get_count(&status, MPI_BYTE, &size);
378  std::vector<char> buffer(size);
379 
380  MPI_Recv(
381  buffer.data(), size, MPI_BYTE, status.MPI_SOURCE, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
382 
383  BufferReader reader(buffer.data());
384  while (reader.getPos() < size) {
385  Interface<D> ifs;
386  reader >> ifs;
387  for (auto& patch : ifs.patches) {
388  auto& ptr = id_to_off_proc_piinfo_map[patch.piinfo->pinfo.id];
389  if (ptr == nullptr) {
390  // need to cast to remove const modifier
391  ptr = std::const_pointer_cast<PatchIfaceInfo<D>>(patch.piinfo);
392  } else {
393  patch.piinfo = ptr;
394  }
395  }
396  rank_id_iface_map.at(rank).at(ifs.id)->merge(ifs);
397  }
398  off_proc_piinfos.clear();
399  for (auto pair : id_to_off_proc_piinfo_map) {
400  off_proc_piinfos.push_back(pair.second);
401  }
402  }
403  // wait for all
404  MPI_Waitall((int)send_requests.size(), &send_requests[0], MPI_STATUSES_IGNORE);
405  }
406 };
407 extern template class Interface<2>;
408 extern template class Interface<3>;
409 } // namespace Schur
410 } // namespace ThunderEgg
411 #endif
ThunderEgg::Schur::Interface::SideTypePiinfo::type
IfaceType< D > type
the IfaceType of this interface on the cooresponding PatchIfaceInfo object.
Definition: Interface.h:78
ThunderEgg::Orthant< D - 1 >
IfaceType.h
IfaceType class.
ThunderEgg::BufferWriter
Class that is used to help serialize objects into a buffer.
Definition: BufferWriter.h:38
ThunderEgg::Schur::Interface::insert
void insert(Side< D > s, std::shared_ptr< const PatchIfaceInfo< D >> piinfo)
Add an associated patch to this interface.
Definition: Interface.h:135
ThunderEgg::Schur::Interface::SideTypePiinfo::getNonConstPiinfo
std::shared_ptr< PatchIfaceInfo< D > > getNonConstPiinfo()
Get the Piinfo object without the const modifer.
Definition: Interface.h:89
ThunderEgg::Schur::Interface::Interface
Interface()=default
Construct a new Interface object.
ThunderEgg::NbrType::Coarse
@ Coarse
The neighbor is at a coarser refinement level.
ThunderEgg::Serializable
Interface for serializing objects.
Definition: Serializable.h:34
ThunderEgg::Schur::Interface::SideTypePiinfo
A struct for each patch associated with this interface.
Definition: Interface.h:68
ThunderEgg::Schur::IfaceType
An enum-style class that represents interface types.
Definition: IfaceType.h:38
ThunderEgg::Schur::Interface::SideTypePiinfo::piinfo
std::shared_ptr< const PatchIfaceInfo< D > > piinfo
the PatchIfaceInfo object
Definition: Interface.h:82
ThunderEgg::Schur::Interface::SideTypePiinfo::SideTypePiinfo
SideTypePiinfo(Side< D > side, IfaceType< D > type, std::shared_ptr< const PatchIfaceInfo< D >> piinfo)
Construct a new Side Type Piinfo object.
Definition: Interface.h:101
ThunderEgg::Schur::Interface::deserialize
int deserialize(char *buffer)
Deserialize an object.
Definition: Interface.h:188
ThunderEgg::Schur::Interface::merge
void merge(Interface< D > ifs)
Add the patches from the given Interface object to this object.
Definition: Interface.h:168
ThunderEgg::NbrType
NbrType
The type of neighbor.
Definition: NbrType.h:34
ThunderEgg
The ThunderEgg namespace.
Definition: BiLinearGhostFiller.h:31
ThunderEgg::BufferReader::getPos
int getPos()
get the current position in the buffer
Definition: BufferReader.h:74
ThunderEgg::Schur::Interface::global_index
int global_index
The global index of the interface.
Definition: Interface.h:61
ThunderEgg::Schur::Interface::patches
std::vector< SideTypePiinfo > patches
SideTypePiinfo structs associated with this interface.
Definition: Interface.h:113
BufferWriter.h
BufferWriter class.
ThunderEgg::NbrType::Fine
@ Fine
The nighbor is at a finer refinement level.
ThunderEgg::NbrType::Normal
@ Normal
The neighbor is at the same refinement level.
ThunderEgg::Schur::Interface::Interface
Interface(int id)
Construct a new Interface object.
Definition: Interface.h:125
ThunderEgg::Schur::Interface::SideTypePiinfo::side
Side< D > side
the side of the PatchIfaceInfo object that this interface is on.
Definition: Interface.h:74
ThunderEgg::Schur::Interface::id
int id
The id of the interface.
Definition: Interface.h:53
ThunderEgg::Schur::Interface::EnumerateIfacesFromPiinfoVector
static void EnumerateIfacesFromPiinfoVector(std::vector< std::shared_ptr< const PatchIfaceInfo< D >>> piinfos, std::map< int, std::map< int, std::shared_ptr< Interface< D >>>> &rank_id_iface_map, std::vector< std::shared_ptr< PatchIfaceInfo< D >>> &off_proc_piinfos)
Will enumerate a map from interface id to this rank's interfaces, will also do any neccesary communic...
Definition: Interface.h:317
ThunderEgg::Schur::Interface::local_index
int local_index
The local index of the interface.
Definition: Interface.h:57
ThunderEgg::BufferWriter::getPos
int getPos()
get the current position in the buffer
Definition: BufferWriter.h:80
ThunderEgg::BufferReader
Class that is used to help read serialized objects from a buffer.
Definition: BufferReader.h:38
BufferReader.h
BufferReader class.
ThunderEgg::Face
Enum-style class for the faces of an n-dimensional cube.
Definition: Face.h:41
ThunderEgg::Schur::Interface::serialize
int serialize(char *buffer) const
Serialize object into buffer.
Definition: Interface.h:174
PatchIfaceInfo.h
PatchIfaceInfo class.
ThunderEgg::RuntimeError
ThunderEgg runtime exception.
Definition: RuntimeError.h:36
ThunderEgg::Schur::Interface
This is a set of PatchIfaceInfo objects for an interface.
Definition: Interface.h:47
ThunderEgg::Schur::PatchIfaceInfo
This decorates a PatchInfo object with a IfaceInfo object for each side of the patch.
Definition: PatchIfaceInfo.h:40