ThunderEgg  1.0.0
PatchIfaceInfo.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_PATCHIFACEINFO_H
22 #define THUNDEREGG_SCHUR_PATCHIFACEINFO_H
23 
32 namespace ThunderEgg {
33 namespace Schur {
39 template<int D>
41 {
42 public:
50  std::array<std::shared_ptr<IfaceInfo<D>>, Side<D>::number_of> iface_info;
54  PatchIfaceInfo() = default;
61  : pinfo(pinfo)
62  {
63  // create iface objects
64  for (Side<D> s : Side<D>::getValues()) {
65  if (pinfo.hasNbr(s)) {
66  switch (pinfo.getNbrType(s)) {
67  case NbrType::Normal:
68  setIfaceInfo(s, std::make_shared<NormalIfaceInfo<D>>(pinfo, s));
69  break;
70  case NbrType::Fine:
71  setIfaceInfo(s, std::make_shared<FineIfaceInfo<D>>(pinfo, s));
72  break;
73  case NbrType::Coarse:
74  setIfaceInfo(s, std::make_shared<CoarseIfaceInfo<D>>(pinfo, s));
75  break;
76  default:
77  throw RuntimeError("Unsupported NbrType");
78  }
79  }
80  }
81  }
88  void setIfaceInfo(Side<D> s, std::shared_ptr<IfaceInfo<D>> info)
89  {
90  iface_info[s.getIndex()] = info;
91  }
98  std::shared_ptr<IfaceInfo<D>> getIfaceInfo(Side<D> s) { return iface_info[s.getIndex()]; }
105  std::shared_ptr<const IfaceInfo<D>> getIfaceInfo(Side<D> s) const
106  {
107  return iface_info[s.getIndex()];
108  }
116  std::shared_ptr<NormalIfaceInfo<D>> getNormalIfaceInfo(Side<D> s)
117  {
118  return std::dynamic_pointer_cast<NormalIfaceInfo<D>>(iface_info[s.getIndex()]);
119  }
127  std::shared_ptr<const NormalIfaceInfo<D>> getNormalIfaceInfo(Side<D> s) const
128  {
129  return std::dynamic_pointer_cast<const NormalIfaceInfo<D>>(iface_info[s.getIndex()]);
130  }
138  std::shared_ptr<CoarseIfaceInfo<D>> getCoarseIfaceInfo(Side<D> s)
139  {
140  return std::dynamic_pointer_cast<CoarseIfaceInfo<D>>(iface_info[s.getIndex()]);
141  }
149  std::shared_ptr<const CoarseIfaceInfo<D>> getCoarseIfaceInfo(Side<D> s) const
150  {
151  return std::dynamic_pointer_cast<const CoarseIfaceInfo<D>>(iface_info[s.getIndex()]);
152  }
160  std::shared_ptr<FineIfaceInfo<D>> getFineIfaceInfo(Side<D> s)
161  {
162  return std::dynamic_pointer_cast<FineIfaceInfo<D>>(iface_info[s.getIndex()]);
163  }
171  std::shared_ptr<const FineIfaceInfo<D>> getFineIfaceInfo(Side<D> s) const
172  {
173  return std::dynamic_pointer_cast<const FineIfaceInfo<D>>(iface_info[s.getIndex()]);
174  }
175  int serialize(char* buffer) const
176  {
177  BufferWriter writer(buffer);
178  writer << pinfo;
179  return writer.getPos();
180  }
181  int deserialize(char* buffer)
182  {
183  BufferReader reader(buffer);
184  auto pinfo_in = std::make_shared<PatchInfo<D>>();
185  reader >> pinfo;
186  // create iface objects
187  for (Side<D> s : Side<D>::getValues()) {
188  if (pinfo.hasNbr(s)) {
189  switch (pinfo.getNbrType(s)) {
190  case NbrType::Normal:
191  setIfaceInfo(s, std::make_shared<NormalIfaceInfo<D>>(pinfo, s));
192  break;
193  case NbrType::Fine:
194  setIfaceInfo(s, std::make_shared<FineIfaceInfo<D>>(pinfo, s));
195  break;
196  case NbrType::Coarse:
197  setIfaceInfo(s, std::make_shared<CoarseIfaceInfo<D>>(pinfo, s));
198  break;
199  default:
200  throw RuntimeError("Unsupported NbrType");
201  }
202  }
203  }
204  return reader.getPos();
205  }
206 };
207 extern template class PatchIfaceInfo<2>;
208 extern template class PatchIfaceInfo<3>;
209 } // namespace Schur
210 } // namespace ThunderEgg
211 #endif
RuntimeError.h
RuntimeError struct.
ThunderEgg::Schur::PatchIfaceInfo::serialize
int serialize(char *buffer) const
Serialize object into buffer.
Definition: PatchIfaceInfo.h:175
ThunderEgg::Schur::PatchIfaceInfo::PatchIfaceInfo
PatchIfaceInfo(const PatchInfo< D > &pinfo)
Construct a new PatchIfaceInfo object.
Definition: PatchIfaceInfo.h:60
ThunderEgg::BufferWriter
Class that is used to help serialize objects into a buffer.
Definition: BufferWriter.h:38
ThunderEgg::Schur::PatchIfaceInfo::getFineIfaceInfo
std::shared_ptr< FineIfaceInfo< D > > getFineIfaceInfo(Side< D > s)
Get the FineIfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:160
ThunderEgg::Schur::PatchIfaceInfo::getFineIfaceInfo
std::shared_ptr< const FineIfaceInfo< D > > getFineIfaceInfo(Side< D > s) const
Get the FineIfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:171
ThunderEgg::NbrType::Coarse
@ Coarse
The neighbor is at a coarser refinement level.
ThunderEgg::Schur::PatchIfaceInfo::getCoarseIfaceInfo
std::shared_ptr< CoarseIfaceInfo< D > > getCoarseIfaceInfo(Side< D > s)
Get the CoarseIfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:138
ThunderEgg::Schur::NormalIfaceInfo
This represents an interface where the neighbor is at the same refinement level.
Definition: NormalIfaceInfo.h:38
ThunderEgg::Serializable
Interface for serializing objects.
Definition: Serializable.h:34
CoarseIfaceInfo.h
CoarseIfaceInfo class.
ThunderEgg::Schur::PatchIfaceInfo::getNormalIfaceInfo
std::shared_ptr< const NormalIfaceInfo< D > > getNormalIfaceInfo(Side< D > s) const
Get the NormalIfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:127
ThunderEgg::Schur::PatchIfaceInfo::getCoarseIfaceInfo
std::shared_ptr< const CoarseIfaceInfo< D > > getCoarseIfaceInfo(Side< D > s) const
Get the CoarseIfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:149
NormalIfaceInfo.h
NormalIfaceInfo class.
ThunderEgg::Schur::CoarseIfaceInfo
Represents the interfaces where the neighbor is at a coarser refinement level.
Definition: CoarseIfaceInfo.h:41
FineIfaceInfo.h
FineIfaceInfo class.
ThunderEgg
The ThunderEgg namespace.
Definition: BiLinearGhostFiller.h:31
ThunderEgg::Face::getIndex
size_t getIndex() const
Get the index for this Face.
Definition: Face.h:452
ThunderEgg::Schur::PatchIfaceInfo::setIfaceInfo
void setIfaceInfo(Side< D > s, std::shared_ptr< IfaceInfo< D >> info)
Set the IfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:88
ThunderEgg::Schur::IfaceInfo
The IfaceInfo class represents the information for an interface on a given side of the patch.
Definition: IfaceInfo.h:43
ThunderEgg::Schur::PatchIfaceInfo::getIfaceInfo
std::shared_ptr< IfaceInfo< D > > getIfaceInfo(Side< D > s)
Get the IfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:98
ThunderEgg::PatchInfo
Contains metadata for a patch.
Definition: PatchInfo.h:51
ThunderEgg::BufferReader::getPos
int getPos()
get the current position in the buffer
Definition: BufferReader.h:74
ThunderEgg::Schur::PatchIfaceInfo::iface_info
std::array< std::shared_ptr< IfaceInfo< D > >, Side< D >::number_of > iface_info
Array of IfaceInfo objects.
Definition: PatchIfaceInfo.h:50
ThunderEgg::Schur::FineIfaceInfo
Represents the interfaces where the neighbors are at a finer refinement level.
Definition: FineIfaceInfo.h:41
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::PatchIfaceInfo::getNormalIfaceInfo
std::shared_ptr< NormalIfaceInfo< D > > getNormalIfaceInfo(Side< D > s)
Get the NormalIfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:116
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::Schur::PatchIfaceInfo::PatchIfaceInfo
PatchIfaceInfo()=default
Construct a new PatchIfaceInfo object.
ThunderEgg::Schur::PatchIfaceInfo::pinfo
PatchInfo< D > pinfo
associated PatchInfo object
Definition: PatchIfaceInfo.h:46
ThunderEgg::Schur::PatchIfaceInfo::getIfaceInfo
std::shared_ptr< const IfaceInfo< D > > getIfaceInfo(Side< D > s) const
Get the IfaceInfo object on a given side of the patch.
Definition: PatchIfaceInfo.h:105
ThunderEgg::RuntimeError
ThunderEgg runtime exception.
Definition: RuntimeError.h:36
ThunderEgg::Schur::PatchIfaceInfo::deserialize
int deserialize(char *buffer)
Deserialize an object.
Definition: PatchIfaceInfo.h:181
ThunderEgg::Schur::PatchIfaceInfo
This decorates a PatchInfo object with a IfaceInfo object for each side of the patch.
Definition: PatchIfaceInfo.h:40