ThunderEgg  1.0.0
PatchInfo.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) 2017-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 #ifndef THUNDEREGG_PATCHINFO_H
21 #define THUNDEREGG_PATCHINFO_H
22 
28 #include <ThunderEgg/FineNbrInfo.h>
30 #include <ThunderEgg/Orthant.h>
33 #include <bitset>
34 
35 namespace ThunderEgg {
50 template<int D>
51 class PatchInfo : public Serializable
52 {
53 private:
58  std::array<std::unique_ptr<NbrInfoBase>, Face<D, D>::sum_of_faces> nbr_infos;
59 
60 public:
65  int id = 0;
69  int local_index = 0;
73  int global_index = 0;
77  int refine_level = -1;
83  int parent_id = -1;
87  int parent_rank = -1;
93  std::array<int, Orthant<D>::num_orthants> child_ids;
99  std::array<int, Orthant<D>::num_orthants> child_ranks;
107  int rank = -1;
117  std::array<int, D> ns;
121  std::array<double, D> starts;
125  std::array<double, D> spacings;
126 
132  {
133  starts.fill(0);
134  ns.fill(1);
135  spacings.fill(1);
136  child_ids.fill(-1);
137  child_ranks.fill(-1);
138  }
139 
145  PatchInfo(const PatchInfo<D>& other_pinfo)
146  : id(other_pinfo.id)
147  , local_index(other_pinfo.local_index)
148  , global_index(other_pinfo.global_index)
149  , refine_level(other_pinfo.refine_level)
150  , parent_id(other_pinfo.parent_id)
151  , parent_rank(other_pinfo.parent_rank)
152  , child_ids(other_pinfo.child_ids)
153  , child_ranks(other_pinfo.child_ranks)
154  , num_ghost_cells(other_pinfo.num_ghost_cells)
155  , rank(other_pinfo.rank)
156  , orth_on_parent(other_pinfo.orth_on_parent)
157  , ns(other_pinfo.ns)
158  , starts(other_pinfo.starts)
159  , spacings(other_pinfo.spacings)
160  {
161  for (size_t i = 0; i < other_pinfo.nbr_infos.size(); i++) {
162  if (other_pinfo.nbr_infos[i] != nullptr) {
163  nbr_infos[i] = other_pinfo.nbr_infos[i]->clone();
164  }
165  }
166  }
173  PatchInfo<D>& operator=(const PatchInfo<D>& other_pinfo)
174  {
175  id = other_pinfo.id;
176  local_index = other_pinfo.local_index;
177  global_index = other_pinfo.global_index;
178  refine_level = other_pinfo.refine_level;
179  parent_id = other_pinfo.parent_id;
180  parent_rank = other_pinfo.parent_rank;
181  child_ids = other_pinfo.child_ids;
182  child_ranks = other_pinfo.child_ranks;
183  num_ghost_cells = other_pinfo.num_ghost_cells;
184  rank = other_pinfo.rank;
185  orth_on_parent = other_pinfo.orth_on_parent;
186  ns = other_pinfo.ns;
187  starts = other_pinfo.starts;
188  spacings = other_pinfo.spacings;
189 
190  for (size_t i = 0; i < other_pinfo.nbr_infos.size(); i++) {
191  if (other_pinfo.nbr_infos[i] != nullptr) {
192  nbr_infos[i] = other_pinfo.nbr_infos[i]->clone();
193  } else {
194  nbr_infos[i] = nullptr;
195  }
196  }
197  return *this;
198  }
207  friend bool operator<(const PatchInfo& l, const PatchInfo& r) { return l.id < r.id; }
208  template<int M>
209  void setNbrInfo(Face<D, M> f, nullptr_t)
210  {
211  nbr_infos[Face<D, M>::sum_of_faces + f.getIndex()] = nullptr;
212  }
220  template<int M>
221  void setNbrInfo(Face<D, M> f, NbrInfo<M>* nbr_info)
222  {
223  nbr_infos[Face<D, M>::sum_of_faces + f.getIndex()].reset(nbr_info);
224  }
231  template<int M>
233  {
234  return nbr_infos[Face<D, M>::sum_of_faces + s.getIndex()]->getNbrType();
235  }
244  template<int M>
246  {
247  return *dynamic_cast<NormalNbrInfo<M>*>(
248  nbr_infos[Face<D, M>::sum_of_faces + s.getIndex()].get());
249  }
257  template<int M>
259  {
260  return *dynamic_cast<CoarseNbrInfo<M>*>(
261  nbr_infos[Face<D, M>::sum_of_faces + s.getIndex()].get());
262  }
271  template<int M>
273  {
274  return *dynamic_cast<FineNbrInfo<M>*>(nbr_infos[Face<D, M>::sum_of_faces + s.getIndex()].get());
275  }
283  template<int M>
284  inline bool hasNbr(Face<D, M> s) const
285  {
286  return nbr_infos[Face<D, M>::sum_of_faces + s.getIndex()] != nullptr;
287  }
291  inline bool hasNbr() const
292  {
293  return std::any_of(
294  nbr_infos.begin(), nbr_infos.end(), [](const std::unique_ptr<NbrInfoBase>& nbr_info) {
295  return nbr_info != nullptr;
296  });
297  }
303  inline bool hasCoarseParent() const { return orth_on_parent != Orthant<D>::null(); }
309  void setNeighborLocalIndexes(const std::map<int, int>& id_to_local_index_map)
310  {
311  for (size_t i = 0; i < nbr_infos.size(); i++) {
312  if (nbr_infos[i] != nullptr) {
313  nbr_infos[i]->setLocalIndexes(id_to_local_index_map);
314  }
315  }
316  }
322  void setNeighborGlobalIndexes(const std::map<int, int>& id_to_global_index_map)
323  {
324  for (size_t i = 0; i < nbr_infos.size(); i++) {
325  if (nbr_infos[i] != nullptr) {
326  nbr_infos[i]->setGlobalIndexes(id_to_global_index_map);
327  }
328  }
329  }
333  std::deque<int> getNbrIds() const
334  {
335  std::deque<int> retval;
336  for (size_t i = 0; i < nbr_infos.size(); i++) {
337  if (nbr_infos[i] != nullptr) {
338  nbr_infos[i]->getNbrIds(retval);
339  }
340  }
341  return retval;
342  }
346  std::deque<int> getNbrRanks() const
347  {
348  std::deque<int> retval;
349  for (size_t i = 0; i < nbr_infos.size(); i++) {
350  if (nbr_infos[i] != nullptr) {
351  nbr_infos[i]->getNbrRanks(retval);
352  }
353  }
354  return retval;
355  }
356  template<int M>
357  void serializeNeighbors(BufferWriter& writer) const
358  {
359  std::bitset<Face<D, M>::number_of> has_nbr;
360  for (Face<D, M> f : Face<D, M>::getValues()) {
361  has_nbr[f.getIndex()] = hasNbr(f);
362  }
363  writer << has_nbr;
364  for (Face<D, M> f : Face<D, M>::getValues()) {
365  if (hasNbr(f)) {
366  NbrType type = getNbrType(f);
367  writer << type;
368  switch (type) {
369  case NbrType::Normal: {
370  NormalNbrInfo<M> tmp = getNormalNbrInfo(f);
371  writer << tmp;
372  } break;
373  case NbrType::Fine: {
374  FineNbrInfo<M> tmp = getFineNbrInfo(f);
375  writer << tmp;
376  } break;
377  case NbrType::Coarse: {
378  CoarseNbrInfo<M> tmp = getCoarseNbrInfo(f);
379  writer << tmp;
380  } break;
381  default:
382  throw RuntimeError("Unsupported NbrType");
383  }
384  }
385  }
386  if constexpr (M > 0) {
387  serializeNeighbors<M - 1>(writer);
388  }
389  }
390  template<int M>
391  void deserializeNeighbors(BufferReader& reader)
392  {
393  std::bitset<Face<D, M>::number_of> has_nbr;
394  reader >> has_nbr;
395  for (Face<D, M> f : Face<D, M>::getValues()) {
396  if (has_nbr[f.getIndex()]) {
397  NbrType type;
398  reader >> type;
399  NbrInfo<M>* info;
400  switch (type) {
401  case NbrType::Normal:
402  info = new NormalNbrInfo<M>();
403  reader >> *static_cast<NormalNbrInfo<M>*>(info);
404  break;
405  case NbrType::Fine:
406  info = new FineNbrInfo<M>();
407  reader >> *static_cast<FineNbrInfo<M>*>(info);
408  break;
409  case NbrType::Coarse:
410  info = new CoarseNbrInfo<M>();
411  reader >> *static_cast<CoarseNbrInfo<M>*>(info);
412  break;
413  default:
414  throw RuntimeError("Unsupported NbrType");
415  }
416 
417  setNbrInfo(f, info);
418  }
419  }
420  if constexpr (M > 0) {
421  deserializeNeighbors<M - 1>(reader);
422  }
423  }
424  int serialize(char* buffer) const
425  {
426  BufferWriter writer(buffer);
427  writer << id;
428  writer << ns;
429  writer << refine_level;
430  writer << parent_id;
431  writer << parent_rank;
432  writer << child_ids;
433  writer << child_ranks;
434  writer << rank;
435  writer << orth_on_parent;
436  writer << starts;
437  writer << spacings;
438 
439  serializeNeighbors<D - 1>(writer);
440 
441  return writer.getPos();
442  }
443  int deserialize(char* buffer)
444  {
445  BufferReader reader(buffer);
446  reader >> id;
447  reader >> ns;
448  reader >> refine_level;
449  reader >> parent_id;
450  reader >> parent_rank;
451  reader >> child_ids;
452  reader >> child_ranks;
453  reader >> rank;
454  reader >> orth_on_parent;
455  reader >> starts;
456  reader >> spacings;
457 
458  deserializeNeighbors<D - 1>(reader);
459 
460  return reader.getPos();
461  }
462 };
463 template<int D>
464 void
465 to_json(tpl::nlohmann::json& j, const PatchInfo<D>& pinfo)
466 {
467  j["id"] = pinfo.id;
468  j["parent_id"] = pinfo.parent_id;
469  j["parent_rank"] = pinfo.parent_rank;
470  j["orth_on_parent"] = pinfo.orth_on_parent;
471  j["rank"] = pinfo.rank;
472  j["starts"] = pinfo.starts;
473  j["lengths"] = pinfo.spacings;
474  j["refine_level"] = pinfo.refine_level;
475  for (int i = 0; i < D; i++) {
476  j["lengths"][i] = pinfo.spacings[i] * pinfo.ns[i];
477  }
478  j["nbrs"] = tpl::nlohmann::json::array();
479  for (Side<D> s : Side<D>::getValues()) {
480  if (pinfo.hasNbr(s)) {
481  switch (pinfo.getNbrType(s)) {
482  case NbrType::Normal:
483  j["nbrs"].push_back(pinfo.getNormalNbrInfo(s));
484  break;
485  case NbrType::Coarse:
486  j["nbrs"].push_back(pinfo.getCoarseNbrInfo(s));
487  break;
488  case NbrType::Fine:
489  j["nbrs"].push_back(pinfo.getFineNbrInfo(s));
490  break;
491  default:
492  throw RuntimeError("Unsupported NbrType");
493  }
494  j["nbrs"].back()["type"] = pinfo.getNbrType(s);
495  j["nbrs"].back()["side"] = s;
496  }
497  }
498  if constexpr (D == 3) {
499  j["edge_nbrs"] = tpl::nlohmann::json::array();
500  for (Edge e : Edge::getValues()) {
501  if (pinfo.hasNbr(e)) {
502  switch (pinfo.getNbrType(e)) {
503  case NbrType::Normal:
504  j["edge_nbrs"].push_back(pinfo.getNormalNbrInfo(e));
505  break;
506  case NbrType::Coarse:
507  j["edge_nbrs"].push_back(pinfo.getCoarseNbrInfo(e));
508  break;
509  case NbrType::Fine:
510  j["edge_nbrs"].push_back(pinfo.getFineNbrInfo(e));
511  break;
512  default:
513  throw RuntimeError("Unsupported NbrType");
514  }
515  j["edge_nbrs"].back()["type"] = pinfo.getNbrType(e);
516  j["edge_nbrs"].back()["edge"] = e;
517  }
518  }
519  }
520  if constexpr (D >= 2) {
521  j["corner_nbrs"] = tpl::nlohmann::json::array();
522  for (Corner<D> c : Corner<D>::getValues()) {
523  if (pinfo.hasNbr(c)) {
524  switch (pinfo.getNbrType(c)) {
525  case NbrType::Normal:
526  j["corner_nbrs"].push_back(pinfo.getNormalNbrInfo(c));
527  break;
528  case NbrType::Coarse:
529  j["corner_nbrs"].push_back(pinfo.getCoarseNbrInfo(c));
530  break;
531  case NbrType::Fine:
532  j["corner_nbrs"].push_back(pinfo.getFineNbrInfo(c));
533  break;
534  default:
535  throw RuntimeError("Unsupported NbrType");
536  }
537  j["corner_nbrs"].back()["type"] = pinfo.getNbrType(c);
538  j["corner_nbrs"].back()["corner"] = c;
539  }
540  }
541  }
542  if (pinfo.child_ids[0] != -1) {
543  j["child_ids"] = pinfo.child_ids;
544  j["child_ranks"] = pinfo.child_ranks;
545  }
546 }
547 template<int D>
548 void
549 from_json(const tpl::nlohmann::json& j, PatchInfo<D>& pinfo)
550 {
551  pinfo.id = j["id"];
552  pinfo.parent_id = j["parent_id"];
553  pinfo.parent_rank = j["parent_rank"];
554  if (j.contains("orth_on_parent")) {
555  j["orth_on_parent"].get_to(pinfo.orth_on_parent);
556  }
557  if (j.contains("refine_level")) {
558  pinfo.refine_level = j["refine_level"];
559  }
560  pinfo.rank = j["rank"];
561  pinfo.starts = j["starts"].get<std::array<double, D>>();
562  pinfo.spacings = j["lengths"].get<std::array<double, D>>();
563  pinfo.ns.fill(1);
564  for (const auto& nbr_j : j["nbrs"]) {
565  Side<D> s = nbr_j["side"].get<Side<D>>();
566  switch (nbr_j["type"].get<NbrType>()) {
567  case NbrType::Normal:
568  pinfo.setNbrInfo(s, new NormalNbrInfo<D - 1>());
569  pinfo.getNormalNbrInfo(s) = nbr_j.get<NormalNbrInfo<D - 1>>();
570  break;
571  case NbrType::Coarse:
572  pinfo.setNbrInfo(s, new CoarseNbrInfo<D - 1>());
573  pinfo.getCoarseNbrInfo(s) = nbr_j.get<CoarseNbrInfo<D - 1>>();
574  break;
575  case NbrType::Fine:
576  pinfo.setNbrInfo(s, new FineNbrInfo<D - 1>());
577  pinfo.getFineNbrInfo(s) = nbr_j.get<FineNbrInfo<D - 1>>();
578  break;
579  default:
580  throw RuntimeError("Unsupported NbrType");
581  }
582  }
583  if constexpr (D == 3) {
584  if (j.contains("edge_nbrs")) {
585  for (const auto& nbr_j : j["edge_nbrs"]) {
586  Edge e = nbr_j["edge"].get<Edge>();
587  switch (nbr_j["type"].get<NbrType>()) {
588  case NbrType::Normal:
589  pinfo.setNbrInfo(e, new NormalNbrInfo<1>());
590  pinfo.getNormalNbrInfo(e) = nbr_j.get<NormalNbrInfo<1>>();
591  break;
592  case NbrType::Coarse:
593  pinfo.setNbrInfo(e, new CoarseNbrInfo<1>());
594  pinfo.getCoarseNbrInfo(e) = nbr_j.get<CoarseNbrInfo<1>>();
595  break;
596  case NbrType::Fine:
597  pinfo.setNbrInfo(e, new FineNbrInfo<1>());
598  pinfo.getFineNbrInfo(e) = nbr_j.get<FineNbrInfo<1>>();
599  break;
600  default:
601  throw RuntimeError("Unsupported NbrType");
602  }
603  }
604  }
605  }
606  if constexpr (D >= 2) {
607  if (j.contains("corner_nbrs")) {
608  for (const auto& nbr_j : j["corner_nbrs"]) {
609  Corner<D> c = nbr_j["corner"].get<Corner<D>>();
610  switch (nbr_j["type"].get<NbrType>()) {
611  case NbrType::Normal:
612  pinfo.setNbrInfo(c, new NormalNbrInfo<0>());
613  pinfo.getNormalNbrInfo(c) = nbr_j.get<NormalNbrInfo<0>>();
614  break;
615  case NbrType::Coarse:
616  pinfo.setNbrInfo(c, new CoarseNbrInfo<0>());
617  pinfo.getCoarseNbrInfo(c) = nbr_j.get<CoarseNbrInfo<0>>();
618  break;
619  case NbrType::Fine:
620  pinfo.setNbrInfo(c, new FineNbrInfo<0>());
621  pinfo.getFineNbrInfo(c) = nbr_j.get<FineNbrInfo<0>>();
622  break;
623  default:
624  throw RuntimeError("Unsupported NbrType");
625  }
626  }
627  }
628  }
629  if (j.contains("child_ids")) {
630  pinfo.child_ids = j["child_ids"].get<std::array<int, Orthant<D>::num_orthants>>();
631  pinfo.child_ranks = j["child_ranks"].get<std::array<int, Orthant<D>::num_orthants>>();
632  }
633 }
634 extern template class PatchInfo<2>;
635 extern template class PatchInfo<3>;
636 } // namespace ThunderEgg
637 #endif
ThunderEgg::PatchInfo::serialize
int serialize(char *buffer) const
Serialize object into buffer.
Definition: PatchInfo.h:424
ThunderEgg::PatchInfo::PatchInfo
PatchInfo()
Construct a new Patch Info object starts, ns, and spacings are all set to 0.
Definition: PatchInfo.h:131
RuntimeError.h
RuntimeError struct.
ThunderEgg::PatchInfo::getFineNbrInfo
FineNbrInfo< M > & getFineNbrInfo(Face< D, M > s) const
Get the FineNbrInfo object.
Definition: PatchInfo.h:272
ThunderEgg::Orthant::null
static Orthant< D > null()
null value
Definition: Orthant.h:68
ThunderEgg::Orthant
An enum-style class that represents the octants of a cube.
Definition: Orthant.h:43
ThunderEgg::PatchInfo::orth_on_parent
Orthant< D > orth_on_parent
The orthant of the parent that this parent resides on.
Definition: PatchInfo.h:113
ThunderEgg::BufferWriter
Class that is used to help serialize objects into a buffer.
Definition: BufferWriter.h:38
ThunderEgg::PatchInfo::getNbrRanks
std::deque< int > getNbrRanks() const
return a vector of neighbor ranks
Definition: PatchInfo.h:346
ThunderEgg::PatchInfo::hasNbr
bool hasNbr(Face< D, M > s) const
Return whether the patch has a neighbor.
Definition: PatchInfo.h:284
ThunderEgg::NbrType::Coarse
@ Coarse
The neighbor is at a coarser refinement level.
ThunderEgg::Serializable
Interface for serializing objects.
Definition: Serializable.h:34
ThunderEgg::PatchInfo::getNormalNbrInfo
NormalNbrInfo< M > & getNormalNbrInfo(Face< D, M > s) const
Get the NormalNbrInfo object for a side.
Definition: PatchInfo.h:245
ThunderEgg::PatchInfo::ns
std::array< int, D > ns
The number of cells in each direction.
Definition: PatchInfo.h:117
ThunderEgg::PatchInfo::parent_id
int parent_id
The id of the parent patch.
Definition: PatchInfo.h:83
ThunderEgg::PatchInfo::operator<
friend bool operator<(const PatchInfo &l, const PatchInfo &r)
Compare the ids of the patches.
Definition: PatchInfo.h:207
ThunderEgg::PatchInfo::setNeighborLocalIndexes
void setNeighborLocalIndexes(const std::map< int, int > &id_to_local_index_map)
Set the local indexes in the NbrInfo objects.
Definition: PatchInfo.h:309
ThunderEgg::PatchInfo::local_index
int local_index
The local index of the patch in the Domain.
Definition: PatchInfo.h:69
ThunderEgg::PatchInfo::operator=
PatchInfo< D > & operator=(const PatchInfo< D > &other_pinfo)
Copy asisgnment.
Definition: PatchInfo.h:173
ThunderEgg::PatchInfo::getNbrType
NbrType getNbrType(Face< D, M > s) const
Get the NbrType for a side.
Definition: PatchInfo.h:232
ThunderEgg::Edge
Face< 3, 1 > Edge
Edge class.
Definition: Face.h:636
ThunderEgg::NbrType
NbrType
The type of neighbor.
Definition: NbrType.h:34
FineNbrInfo.h
FineNbrInfo class.
ThunderEgg::PatchInfo::PatchInfo
PatchInfo(const PatchInfo< D > &other_pinfo)
Copy constructor.
Definition: PatchInfo.h:145
ThunderEgg
The ThunderEgg namespace.
Definition: BiLinearGhostFiller.h:31
ThunderEgg::PatchInfo::id
int id
The globally unique ID of the patch This ID only needs to be unique within a Domain.
Definition: PatchInfo.h:65
ThunderEgg::NbrInfo
Represents information about a patch's neighbor.
Definition: NbrInfo.h:36
ThunderEgg::Face::getIndex
size_t getIndex() const
Get the index for this Face.
Definition: Face.h:452
ThunderEgg::CoarseNbrInfo
Represents a neighbor that is at a coarser refinement level.
Definition: CoarseNbrInfo.h:39
NormalNbrInfo.h
NormalNbrType class.
ThunderEgg::PatchInfo::setNeighborGlobalIndexes
void setNeighborGlobalIndexes(const std::map< int, int > &id_to_global_index_map)
Set the global indexes in the NbrInfo objects.
Definition: PatchInfo.h:322
ThunderEgg::PatchInfo
Contains metadata for a patch.
Definition: PatchInfo.h:51
ThunderEgg::PatchInfo::hasNbr
bool hasNbr() const
Return if this patch has a neighbor.
Definition: PatchInfo.h:291
Serializable.h
Serializable class.
ThunderEgg::PatchInfo::deserialize
int deserialize(char *buffer)
Deserialize an object.
Definition: PatchInfo.h:443
ThunderEgg::PatchInfo::getNbrIds
std::deque< int > getNbrIds() const
return a vector of neighbor ids
Definition: PatchInfo.h:333
ThunderEgg::PatchInfo::rank
int rank
MPI rank of this patch.
Definition: PatchInfo.h:107
ThunderEgg::PatchInfo::setNbrInfo
void setNbrInfo(Face< D, M > f, NbrInfo< M > *nbr_info)
Set the Nbr Info object on a given face.
Definition: PatchInfo.h:221
ThunderEgg::BufferReader::getPos
int getPos()
get the current position in the buffer
Definition: BufferReader.h:74
CoarseNbrInfo.h
CoarseNbrInfo 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::PatchInfo::spacings
std::array< double, D > spacings
The cell spacings in each direction.
Definition: PatchInfo.h:125
ThunderEgg::PatchInfo::child_ids
std::array< int, Orthant< D >::num_orthants > child_ids
The id's of the children.
Definition: PatchInfo.h:93
ThunderEgg::PatchInfo::starts
std::array< double, D > starts
The lower-left-bottom index of the patch.
Definition: PatchInfo.h:121
ThunderEgg::PatchInfo::refine_level
int refine_level
The refinement level.
Definition: PatchInfo.h:77
ThunderEgg::PatchInfo::num_ghost_cells
int num_ghost_cells
Number of ghost cells on each side of the patch.
Definition: PatchInfo.h:103
ThunderEgg::PatchInfo::global_index
int global_index
The global index of the patch in the Domain.
Definition: PatchInfo.h:73
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
ThunderEgg::Face
Enum-style class for the faces of an n-dimensional cube.
Definition: Face.h:41
ThunderEgg::PatchInfo::child_ranks
std::array< int, Orthant< D >::num_orthants > child_ranks
The ranks of the children.
Definition: PatchInfo.h:99
ThunderEgg::FineNbrInfo
Represents neighbors that are at a finer refinement level.
Definition: FineNbrInfo.h:39
ThunderEgg::NormalNbrInfo
Represents a neighbor that is at the same refinement level.
Definition: NormalNbrInfo.h:38
ThunderEgg::PatchInfo::parent_rank
int parent_rank
the rank that the parent patch resides on
Definition: PatchInfo.h:87
ThunderEgg::PatchInfo::getCoarseNbrInfo
CoarseNbrInfo< M > & getCoarseNbrInfo(Face< D, M > s) const
Get the CoarseNbrInfo object.
Definition: PatchInfo.h:258
ThunderEgg::PatchInfo::hasCoarseParent
bool hasCoarseParent() const
Return whether the patch has a coarser parent.
Definition: PatchInfo.h:303
Orthant.h
Orthant class.
ThunderEgg::Face::getValues
static Range getValues()
Get a range of values that can be iterated over.
Definition: Face.h:445