v0.15.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 SETERRQ(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 reference 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 = dynamic_cast<IFACE *>(ptr)))
98 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
99 "Cast impossible %s -> %s",
100 (boost::typeindex::type_id_runtime(*ptr).pretty_name()).c_str(),
101 (boost::typeindex::type_id<IFACE>().pretty_name()).c_str());
103 }
104
105 /**
106 * @brief Get interface pointer to pointer of interface
107 *
108 * \code
109 * // Create moab database
110 * moab::Core mb_instance;
111 * // Access moab database by interface
112 * moab::Interface &moab = mb_instance;
113 *
114 * // Create MoFEM database
115 * MoFEM::Core core(moab);
116 * // Acces MoFEM database by Interface
117 * MoFEM::Interface &m_field = core;
118 *
119 * // Get interface
120 * // Get simple interface is simplified version enabling quick and
121 * // easy construction of problem.
122 * Simple *simple_interface;
123 * // Query interface and get pointer to Simple interface
124 * CHKERR m_field.getInterface(&simple_interface);
125 *
126 * \endcode
127 *
128 *
129 * @param iface const pointer to a interface pointer
130 * @return MoFEMErrorCode
131 */
132 template <class IFACE>
133 inline MoFEMErrorCode getInterface(IFACE **const iface) const {
134 return getInterface<IFACE>(*iface);
135 }
136
137 /**
138 * @brief Get interface pointer to pointer of interface
139 *
140 * \code
141 * // Create moab database
142 * moab::Core mb_instance;
143 * // Access moab database by interface
144 * moab::Interface &moab = mb_instance;
145 *
146 * // Create MoFEM database
147 * MoFEM::Core core(moab);
148 * // Acces MoFEM database by Interface
149 * MoFEM::Interface &m_field = core;
150 *
151 * // Get interface
152 * // Get simple interface is simplified version enabling quick and
153 * // easy construction of problem.
154 * Simple *simple_interface = m_field.getInterface<Simple*>();
155 *
156 * \endcode
157 *
158 * @return IFACE*
159 */
160 template <class IFACE,
161 typename boost::enable_if<boost::is_pointer<IFACE>, int>::type = 0>
162 inline IFACE getInterface() const {
163 typedef typename boost::remove_pointer<IFACE>::type IFaceType;
164 IFaceType *iface = NULL;
165 ierr = getInterface<IFACE>(iface);
166 CHKERRABORT(PETSC_COMM_SELF, ierr);
167 return iface;
168 }
169
170 /**
171 * @brief Get reference to interface
172 *
173 * \code
174 * // Create moab database
175 * moab::Core mb_instance;
176 * // Access moab database by interface
177 * moab::Interface &moab = mb_instance;
178 *
179 * // Create MoFEM database
180 * MoFEM::Core core(moab);
181 * // Acces MoFEM database by Interface
182 * MoFEM::Interface &m_field = core;
183 *
184 * // Get interface
185 * // Get simple interface is simplified version enabling quick and
186 * // easy construction of problem.
187 * Simple &simple_interface = m_field.getInterface<Simple&>();
188 *
189 * \endcode
190 *
191 * @return IFACE&
192 */
193 template <class IFACE, typename boost::enable_if<boost::is_reference<IFACE>,
194 int>::type = 0>
195 inline IFACE getInterface() const {
196 typedef typename boost::remove_reference<IFACE>::type IFaceType;
197 IFaceType *iface = NULL;
198 ierr = getInterface<IFaceType>(iface);
199 CHKERRABORT(PETSC_COMM_SELF, ierr);
200 return *iface;
201 }
202
203 /**
204 * @brief Function returning pointer to interface
205 *
206 * \code
207 * // Create moab database
208 * moab::Core mb_instance;
209 * // Access moab database by interface
210 * moab::Interface &moab = mb_instance;
211 *
212 * // Create MoFEM database
213 * MoFEM::Core core(moab);
214 * // Acces MoFEM database by Interface
215 * MoFEM::Interface &m_field = core;
216 *
217 * // Get interface
218 * // Get simple interface is simplified version enabling quick and
219 * // easy construction of problem.
220 * Simple *simple_interface = m_field.getInterface<Simple>();
221 *
222 * \endcode
223 *
224 * @return IFACE*
225 */
226 template <class IFACE> inline IFACE *getInterface() const {
227 IFACE *iface = NULL;
228 ierr = getInterface<IFACE>(iface);
229 CHKERRABORT(PETSC_COMM_SELF, ierr);
230 return iface;
231 }
232
233 virtual ~UnknownInterface() = default;
234
235 /**
236 * \brief Get library version
237 *
238 * This is library version.
239 *
240 * @return error code
241 */
242 static MoFEMErrorCode getLibVersion(Version &version);
243
244 /**
245 * \brief Get database major version
246 *
247 * This is database version. MoFEM can read DataBase from file created by
248 * older version. Then library version and database version could be
249 * different.
250 *
251 * @return error code
252 */
253 static MoFEMErrorCode getFileVersion(moab::Interface &moab, Version &version);
254
255 /**
256 * \brief Get database major version
257 *
258 * This is database version. MoFEM can read DataBase from file created by
259 * older version. Then library version and database version could be
260 * different.
261 *
262 * @return error code
263 */
265 moab::Interface &moab,
266 Version version = Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR,
267 MoFEM_VERSION_BUILD));
268
269 /**
270 * \brief Get database major version
271 *
272 * Implementation of particular interface could be different than main lib.
273 * For example user could use older interface, to keep back compatibility.
274 *
275 * @return error code
276 */
278
279protected:
280 struct NotKnownClass {};
281
282private:
283 struct UIdTypeMap {
284 boost::typeindex::type_index classIdx;
285 UIdTypeMap(const boost::typeindex::type_index &idx) : classIdx(idx) {}
286 };
287
288 /// Data structure for interfaces Id and class names
289 typedef multi_index_container<
290 UIdTypeMap,
291 indexed_by<
292
293 hashed_unique<member<UIdTypeMap, boost::typeindex::type_index,
295
296 >>
298
300 iFaceTypeMap; ///< Maps implementation to interface type name
301};
302
303} // namespace MoFEM
304
305#endif // __MOFEMUNKNOWNINTERFACE_HPP__
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
const double v
phase velocity of light in medium (cm/ns)
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
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 reference 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)