v0.14.0
Loading...
Searching...
No Matches
UnknownInterface.hpp
Go to the documentation of this file.
1/** \file UnknownInterface.hpp
2 * \brief MoFEM interface
3 *
4 * Low level data structures not used directly by user
5 */
6
7#ifndef __MOFEMUNKNOWNINTERFACE_HPP__
8#define __MOFEMUNKNOWNINTERFACE_HPP__
9
10namespace MoFEM {
11
12struct Version {
17 : majorVersion(MoFEM_VERSION_MAJOR), minorVersion(MoFEM_VERSION_MINOR),
18 buildVersion(MoFEM_VERSION_BUILD) {}
19 Version(const int v[3])
20 : majorVersion(v[0]), minorVersion(v[1]), buildVersion(v[2]) {}
21 Version(const int minor, const int major, const int build)
22 : majorVersion(minor), minorVersion(major), buildVersion(build) {}
23
24 std::string strVersion() {
25 auto str = [](auto v) { return boost::lexical_cast<std::string>(v); };
26 return str(majorVersion) + "." + str(minorVersion) + "." +
27 str(buildVersion);
28 }
29};
30
31/** \brief base class for all interface classes
32 * \ingroup mofem
33 **/
35
36 virtual MoFEMErrorCode
37 query_interface(boost::typeindex::type_index type_index,
38 UnknownInterface **iface) const = 0;
39
40 /**
41 * @brief Register interface
42 *
43 * Example:
44 * \code
45 * ierr = regSubInterface<Simple>(IDD_MOFEMSimple);
46 * CHKERRABORT(PETSC_COMM_SELF, ierr);
47 * \endcode
48 *
49 * @param uuid
50 * @param true
51 * @return MoFEMErrorCode
52 */
53 template <class IFACE>
54 MoFEMErrorCode registerInterface(bool error_if_registration_failed = true) {
56 auto p =
57 iFaceTypeMap.insert(UIdTypeMap(boost::typeindex::type_id<IFACE>()));
58 if (error_if_registration_failed && (!p.second)) {
59 SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
60 "Registration of interface typeid(IFACE).name() = %s failed",
61 typeid(IFACE).name());
62 }
64 }
65
66 /**
67 * @brief Get interface refernce to pointer of interface
68 *
69 * \code
70 * // Create moab database
71 * moab::Core mb_instance;
72 * // Access moab database by interface
73 * moab::Interface &moab = mb_instance;
74 *
75 * // Create MoFEM database
76 * MoFEM::Core core(moab);
77 * // Acces MoFEM database by Interface
78 * MoFEM::Interface &m_field = core;
79 *
80 * // Get interface
81 * // Get simple interface is simplified version enabling quick and
82 * // easy construction of problem.
83 * Simple *simple_interface;
84 * // Query interface and get pointer to Simple interface
85 * CHKERR m_field.getInterface(simple_interface);
86 *
87 * \endcode
88 *
89 * @param iface reference to a interface pointer
90 * @return MoFEMErrorCode
91 */
92 template <class IFACE>
93 inline MoFEMErrorCode getInterface(IFACE *&iface) const {
96 CHKERR query_interface(boost::typeindex::type_id<IFACE>(), &ptr);
97 if (!(iface = static_cast<IFACE *>(ptr)))
98 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "Cast Impossible");
100 }
101
102 /**
103 * @brief Get interface pointer to pointer of interface
104 *
105 * \code
106 * // Create moab database
107 * moab::Core mb_instance;
108 * // Access moab database by interface
109 * moab::Interface &moab = mb_instance;
110 *
111 * // Create MoFEM database
112 * MoFEM::Core core(moab);
113 * // Acces MoFEM database by Interface
114 * MoFEM::Interface &m_field = core;
115 *
116 * // Get interface
117 * // Get simple interface is simplified version enabling quick and
118 * // easy construction of problem.
119 * Simple *simple_interface;
120 * // Query interface and get pointer to Simple interface
121 * CHKERR m_field.getInterface(&simple_interface);
122 *
123 * \endcode
124 *
125 *
126 * @param iface const pointer to a interface pointer
127 * @return MoFEMErrorCode
128 */
129 template <class IFACE>
130 inline MoFEMErrorCode getInterface(IFACE **const iface) const {
131 return getInterface<IFACE>(*iface);
132 }
133
134 /**
135 * @brief Get interface pointer to pointer of interface
136 *
137 * \code
138 * // Create moab database
139 * moab::Core mb_instance;
140 * // Access moab database by interface
141 * moab::Interface &moab = mb_instance;
142 *
143 * // Create MoFEM database
144 * MoFEM::Core core(moab);
145 * // Acces MoFEM database by Interface
146 * MoFEM::Interface &m_field = core;
147 *
148 * // Get interface
149 * // Get simple interface is simplified version enabling quick and
150 * // easy construction of problem.
151 * Simple *simple_interface = m_field.getInterface<Simple*>();
152 *
153 * \endcode
154 *
155 * @return IFACE*
156 */
157 template <class IFACE,
158 typename boost::enable_if<boost::is_pointer<IFACE>, int>::type = 0>
159 inline IFACE getInterface() const {
160 typedef typename boost::remove_pointer<IFACE>::type IFaceType;
161 IFaceType *iface = NULL;
162 ierr = getInterface<IFACE>(iface);
163 CHKERRABORT(PETSC_COMM_SELF, ierr);
164 return iface;
165 }
166
167 /**
168 * @brief Get reference to interface
169 *
170 * \code
171 * // Create moab database
172 * moab::Core mb_instance;
173 * // Access moab database by interface
174 * moab::Interface &moab = mb_instance;
175 *
176 * // Create MoFEM database
177 * MoFEM::Core core(moab);
178 * // Acces MoFEM database by Interface
179 * MoFEM::Interface &m_field = core;
180 *
181 * // Get interface
182 * // Get simple interface is simplified version enabling quick and
183 * // easy construction of problem.
184 * Simple &simple_interface = m_field.getInterface<Simple&>();
185 *
186 * \endcode
187 *
188 * @return IFACE&
189 */
190 template <class IFACE, typename boost::enable_if<boost::is_reference<IFACE>,
191 int>::type = 0>
192 inline IFACE getInterface() const {
193 typedef typename boost::remove_reference<IFACE>::type IFaceType;
194 IFaceType *iface = NULL;
195 ierr = getInterface<IFaceType>(iface);
196 CHKERRABORT(PETSC_COMM_SELF, ierr);
197 return *iface;
198 }
199
200 /**
201 * @brief Function returning pointer to interface
202 *
203 * \code
204 * // Create moab database
205 * moab::Core mb_instance;
206 * // Access moab database by interface
207 * moab::Interface &moab = mb_instance;
208 *
209 * // Create MoFEM database
210 * MoFEM::Core core(moab);
211 * // Acces MoFEM database by Interface
212 * MoFEM::Interface &m_field = core;
213 *
214 * // Get interface
215 * // Get simple interface is simplified version enabling quick and
216 * // easy construction of problem.
217 * Simple *simple_interface = m_field.getInterface<Simple>();
218 *
219 * \endcode
220 *
221 * @return IFACE*
222 */
223 template <class IFACE> inline IFACE *getInterface() const {
224 IFACE *iface = NULL;
225 ierr = getInterface<IFACE>(iface);
226 CHKERRABORT(PETSC_COMM_SELF, ierr);
227 return iface;
228 }
229
230 virtual ~UnknownInterface() = default;
231
232 /**
233 * \brief Get library version
234 *
235 * This is library version.
236 *
237 * @return error code
238 */
239 static MoFEMErrorCode getLibVersion(Version &version);
240
241 /**
242 * \brief Get database major version
243 *
244 * This is database version. MoFEM can read DataBase from file created by
245 * older version. Then library version and database version could be
246 * different.
247 *
248 * @return error code
249 */
250 static MoFEMErrorCode getFileVersion(moab::Interface &moab, Version &version);
251
252 /**
253 * \brief Get database major version
254 *
255 * This is database version. MoFEM can read DataBase from file created by
256 * older version. Then library version and database version could be
257 * different.
258 *
259 * @return error code
260 */
262 moab::Interface &moab,
263 Version version = Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR,
264 MoFEM_VERSION_BUILD));
265
266 /**
267 * \brief Get database major version
268 *
269 * Implementation of particular interface could be different than main lib.
270 * For example user could use older interface, to keep back compatibility.
271 *
272 * @return error code
273 */
275
276protected:
277 struct NotKnownClass {};
278
279private:
280 struct UIdTypeMap {
281 boost::typeindex::type_index classIdx;
282 UIdTypeMap(const boost::typeindex::type_index &idx) : classIdx(idx) {}
283 };
284
285 /// Data structure for interfaces Id and class names
286 typedef multi_index_container<
287 UIdTypeMap,
288 indexed_by<
289
290 hashed_unique<member<UIdTypeMap, boost::typeindex::type_index,
292
293 >>
295
297 iFaceTypeMap; ///< Maps implementation to interface type name
298};
299
300} // namespace MoFEM
301
302#endif // __MOFEMUNKNOWNINTERFACE_HPP__
static Index< 'p', 3 > p
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
const double v
phase velocity of light in medium (cm/ns)
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
boost::typeindex::type_index classIdx
UIdTypeMap(const boost::typeindex::type_index &idx)
base class for all interface classes
multi_index_container< UIdTypeMap, indexed_by< hashed_unique< member< UIdTypeMap, boost::typeindex::type_index, &UIdTypeMap::classIdx > > > > iFaceTypeMap_multiIndex
Data structure for interfaces Id and class names.
iFaceTypeMap_multiIndex iFaceTypeMap
Maps implementation to interface type name.
IFACE * getInterface() const
Function returning pointer to interface.
virtual MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
static MoFEMErrorCode getLibVersion(Version &version)
Get library version.
MoFEMErrorCode registerInterface(bool error_if_registration_failed=true)
Register interface.
MoFEMErrorCode getInterface(IFACE **const iface) const
Get interface pointer to pointer of interface.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
static MoFEMErrorCode setFileVersion(moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
Get database major version.
IFACE getInterface() const
Get interface pointer to pointer of interface.
static MoFEMErrorCode getInterfaceVersion(Version &version)
Get database major version.
static MoFEMErrorCode getFileVersion(moab::Interface &moab, Version &version)
Get database major version.
virtual ~UnknownInterface()=default
std::string strVersion()
Version(const int v[3])
Version(const int minor, const int major, const int build)