Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 :
4 : // This library is free software; you can redistribute it and/or
5 : // modify it under the terms of the GNU Lesser General Public
6 : // License as published by the Free Software Foundation; either
7 : // version 2.1 of the License, or (at your option) any later version.
8 :
9 : // This library is distributed in the hope that it will be useful,
10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 : // Lesser General Public License for more details.
13 :
14 : // You should have received a copy of the GNU Lesser General Public
15 : // License along with this library; if not, write to the Free Software
16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 :
18 :
19 :
20 : #ifndef LIBMESH_LIBMESH_COMMON_H
21 : #define LIBMESH_LIBMESH_COMMON_H
22 :
23 : // These flags should never be used together. -DDEBUG means "turn on
24 : // all the ridiculously expensive error checking"; -DNDEBUG means
25 : // "turn off all the somewhat affordable error checking"
26 : #if defined(DEBUG) && defined(NDEBUG)
27 : # error DEBUG and NDEBUG should never be defined simultaneously
28 : #endif
29 :
30 : // The library configuration options
31 : #include "libmesh/libmesh_config.h"
32 :
33 : // Use actual timestamps or constant dummies (to aid ccache)
34 : #ifdef LIBMESH_ENABLE_TIMESTAMPS
35 : # define LIBMESH_TIME __TIME__
36 : # define LIBMESH_DATE __DATE__
37 : #else
38 : # define LIBMESH_TIME "notime"
39 : # define LIBMESH_DATE "nodate"
40 : #endif
41 :
42 : // C/C++ includes everyone should know about
43 : #include <cstdlib>
44 : #ifdef __PGI
45 : // BSK, Thu Feb 20 08:32:06 CST 2014 - For some reason, unless PGI gets
46 : // <cmath> early this nonsense shows up:
47 : // "/software/x86_64/pgi/12.9/linux86-64/12.9/include/CC/cmath", line 57: error:
48 : // the global scope has no "abs"
49 : // using _STLP_VENDOR_CSTD::abs;
50 : // So include <cmath> as early as possible under the PGI compilers.
51 : # include <cmath>
52 : #endif
53 : #include <complex>
54 : #include <typeinfo> // std::bad_cast
55 : #include <type_traits> // std::decay
56 : #include <functional> // std::less, etc
57 : #include <string>
58 :
59 : // Include the MPI definition
60 : #ifdef LIBMESH_HAVE_MPI
61 : # include "libmesh/ignore_warnings.h"
62 : # include <mpi.h>
63 : # include "libmesh/restore_warnings.h"
64 : #endif
65 :
66 : // Quad precision if we need it
67 : #ifdef LIBMESH_DEFAULT_QUADRUPLE_PRECISION
68 : #include "libmesh/float128_shims.h"
69 : #endif
70 :
71 : // _basic_ library functionality
72 : #include "libmesh/libmesh_base.h"
73 : #include "libmesh/libmesh_exceptions.h"
74 :
75 : // Proxy class for libMesh::out/err output
76 : #include "libmesh/ostream_proxy.h"
77 :
78 : // Make sure the libmesh_nullptr define is available for backwards
79 : // compatibility, although we no longer use it in the library.
80 : #include "libmesh/libmesh_nullptr.h"
81 :
82 : // C++ headers
83 : #include <iomanip> // setprecision, in assertion macros
84 :
85 : namespace libMesh
86 : {
87 : namespace Threads
88 : {
89 : // For thread-safe error-messaging. Definitions in threads.h
90 : void lock_singleton_spin_mutex();
91 : void unlock_singleton_spin_mutex();
92 : }
93 :
94 : // Let's define a couple output streams - these will default
95 : // to cout/cerr, but LibMeshInit (or the user) can also set them to
96 : // something more sophisticated.
97 : //
98 : // We use a proxy class rather than references so they can be
99 : // reseated at runtime.
100 :
101 : extern OStreamProxy out;
102 : extern OStreamProxy err;
103 :
104 : // A namespace for functions used in the bodies of the macros below.
105 : // The macros generally call these functions with __FILE__, __LINE__,
106 : // __DATE__, and __TIME__ in the appropriate order. These should not
107 : // be called by users directly! The implementations can be found in
108 : // libmesh_common.C.
109 : namespace MacroFunctions
110 : {
111 : void here(const char * file, int line, const char * date, const char * time, std::ostream & os = libMesh::err);
112 : void stop(const char * file, int line, const char * date, const char * time);
113 : void report_error(const char * file, int line, const char * date, const char * time, std::ostream & os = libMesh::err);
114 : }
115 :
116 : // Undefine any existing macros
117 : #ifdef Real
118 : # undef Real
119 : #endif
120 :
121 : //#ifdef REAL
122 : //# undef REAL
123 : //#endif
124 :
125 : #ifdef Complex
126 : # undef Complex
127 : #endif
128 :
129 : #ifdef COMPLEX
130 : # undef COMPLEX
131 : #endif
132 :
133 : // Check to see if TOLERANCE has been defined by another
134 : // package, if so we might want to change the name...
135 : #ifdef TOLERANCE
136 : DIE A HORRIBLE DEATH HERE...
137 : # undef TOLERANCE
138 : #endif
139 :
140 :
141 :
142 : // Define the type to use for real numbers
143 :
144 : typedef LIBMESH_DEFAULT_SCALAR_TYPE Real;
145 :
146 : // Define a corresponding tolerance. This is what should be
147 : // considered "good enough" when doing floating point comparisons.
148 : // For example, v == 0 is changed to std::abs(v) < TOLERANCE.
149 :
150 : #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
151 : static constexpr Real TOLERANCE = 2.5e-3;
152 : # if defined (LIBMESH_DEFAULT_TRIPLE_PRECISION) || \
153 : defined (LIBMESH_DEFAULT_QUADRUPLE_PRECISION)
154 : # error Cannot define multiple precision levels
155 : # endif
156 : #endif
157 :
158 : #ifdef LIBMESH_DEFAULT_TRIPLE_PRECISION
159 : static constexpr Real TOLERANCE = 1.e-8;
160 : # if defined (LIBMESH_DEFAULT_QUADRUPLE_PRECISION)
161 : # error Cannot define multiple precision levels
162 : # endif
163 : #endif
164 :
165 : #ifdef LIBMESH_DEFAULT_QUADRUPLE_PRECISION
166 : static constexpr Real TOLERANCE = 1.e-11;
167 : #endif
168 :
169 : #if !defined (LIBMESH_DEFAULT_SINGLE_PRECISION) && \
170 : !defined (LIBMESH_DEFAULT_TRIPLE_PRECISION) && \
171 : !defined (LIBMESH_DEFAULT_QUADRUPLE_PRECISION)
172 : static constexpr Real TOLERANCE = 1.e-6;
173 : #endif
174 :
175 : // Define the type to use for complex numbers
176 : // Always use std::complex<double>, as required by Petsc?
177 : // If your version of Petsc doesn't support
178 : // std::complex<other_precision>, then you'd better just leave
179 : // Real==double
180 : typedef std::complex<Real> Complex;
181 : typedef std::complex<Real> COMPLEX;
182 :
183 :
184 : // Helper functions for complex/real numbers
185 : // to clean up #ifdef LIBMESH_USE_COMPLEX_NUMBERS elsewhere
186 20240710 : template<typename T> inline T libmesh_real(T a) { return a; }
187 : template<typename T> inline T libmesh_imag(T /*a*/) { return 0; }
188 36672082 : template<typename T> inline T libmesh_conj(T a) { return a; }
189 :
190 : template<typename T>
191 0 : inline T libmesh_real(std::complex<T> a) { return std::real(a); }
192 :
193 : template<typename T>
194 : inline T libmesh_imag(std::complex<T> a) { return std::imag(a); }
195 :
196 : template<typename T>
197 : inline std::complex<T> libmesh_conj(std::complex<T> a) { return std::conj(a); }
198 :
199 : // std::isnan() is in <cmath> as of C++11.
200 : template <typename T>
201 111205 : inline bool libmesh_isnan(T x) { return std::isnan(x); }
202 :
203 : template <typename T>
204 : inline bool libmesh_isnan(std::complex<T> a)
205 42952 : { return (std::isnan(std::real(a)) || std::isnan(std::imag(a))); }
206 :
207 : // std::isinf() is in <cmath> as of C++11.
208 : template <typename T>
209 650 : inline bool libmesh_isinf(T x) { return std::isinf(x); }
210 :
211 : template <typename T>
212 : inline bool libmesh_isinf(std::complex<T> a)
213 : { return (std::isinf(std::real(a)) || std::isinf(std::imag(a))); }
214 :
215 : // Define the value type for unknowns in simulations.
216 : // This is either Real or Complex, depending on how
217 : // the library was configures
218 : #if defined (LIBMESH_USE_REAL_NUMBERS)
219 : typedef Real Number;
220 : #elif defined (LIBMESH_USE_COMPLEX_NUMBERS)
221 : typedef Complex Number;
222 : #else
223 : DIE A HORRIBLE DEATH HERE...
224 : #endif
225 :
226 :
227 : // Define the value type for error estimates.
228 : // Since AMR/C decisions don't have to be precise,
229 : // we default to float for memory efficiency.
230 : typedef float ErrorVectorReal;
231 : #define MPI_ERRORVECTORREAL MPI_FLOAT
232 :
233 :
234 : #ifdef LIBMESH_HAVE_MPI
235 :
236 : /**
237 : * MPI Communicator used to initialize libMesh.
238 : */
239 : extern MPI_Comm GLOBAL_COMM_WORLD;
240 : #else
241 :
242 : /**
243 : * Something to use with CHKERRABORT if we're just using PETSc's MPI
244 : * "uni" stub.
245 : */
246 : extern int GLOBAL_COMM_WORLD;
247 : #endif
248 :
249 : // This global variable is to help us deprecate AutoPtr. We can't
250 : // just use libmesh_deprecated() because then you get one print out
251 : // per template instantiation, instead of one total print out.
252 : extern bool warned_about_auto_ptr;
253 :
254 : // These are useful macros that behave like functions in the code.
255 : // If you want to make sure you are accessing a section of code just
256 : // stick a libmesh_here(); in it, for example
257 : #define libmesh_here() \
258 : do { \
259 : libMesh::MacroFunctions::here(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
260 : } while (0)
261 :
262 : // the libmesh_stop() macro will stop the code until a SIGCONT signal
263 : // is received. This is useful, for example, when determining the
264 : // memory used by a given operation. A libmesh_stop() could be
265 : // inserted before and after a questionable operation and the delta
266 : // memory can be obtained from a ps or top. This macro only works for
267 : // serial cases.
268 : #define libmesh_stop() \
269 : do { \
270 : libMesh::MacroFunctions::stop(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
271 : } while (0)
272 :
273 : // The libmesh_dbg_var() macro indicates that an argument to a function
274 : // is used only in debug mode (i.e., when NDEBUG is not defined).
275 : #ifndef NDEBUG
276 : #define libmesh_dbg_var(var) var
277 : #else
278 : #define libmesh_dbg_var(var)
279 : #endif
280 :
281 : // The libmesh_inf_var() macro indicates that an argument to a function
282 : // is used only when infinite elements are enabled
283 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
284 : #define libmesh_inf_var(var) var
285 : #else
286 : #define libmesh_inf_var(var)
287 : #endif
288 :
289 : // The libmesh_assert() macro acts like C's assert(), but throws a
290 : // libmesh_error() (including stack trace, etc) instead of just exiting
291 : #ifdef NDEBUG
292 :
293 : #define libmesh_assert_msg(asserted, msg) ((void) 0)
294 : #define libmesh_exceptionless_assert_msg(asserted, msg) ((void) 0)
295 : #define libmesh_assert_equal_to_msg(expr1,expr2, msg) ((void) 0)
296 : #define libmesh_assert_not_equal_to_msg(expr1,expr2, msg) ((void) 0)
297 : #define libmesh_assert_less_msg(expr1,expr2, msg) ((void) 0)
298 : #define libmesh_assert_greater_msg(expr1,expr2, msg) ((void) 0)
299 : #define libmesh_assert_less_equal_msg(expr1,expr2, msg) ((void) 0)
300 : #define libmesh_assert_greater_equal_msg(expr1,expr2, msg) ((void) 0)
301 :
302 : #else
303 :
304 : #define libmesh_assertion_types(expr1,expr2) \
305 : typedef typename std::decay<decltype(expr1)>::type libmesh_type1; \
306 : typedef typename std::decay<decltype(expr2)>::type libmesh_type2
307 :
308 : #define libmesh_assert_msg(asserted, msg) \
309 : do { \
310 : if (!(asserted)) { \
311 : libmesh_error_msg(msg); \
312 : } } while (0)
313 :
314 : #define libmesh_exceptionless_assert_msg(asserted, msg) \
315 : do { \
316 : if (!(asserted)) { \
317 : libMesh::Threads::lock_singleton_spin_mutex(); \
318 : libMesh::err << "Assertion `" #asserted "' failed." << std::endl; \
319 : libMesh::Threads::unlock_singleton_spin_mutex(); \
320 : libmesh_exceptionless_error(); \
321 : } } while (0)
322 :
323 : #define libmesh_assert_equal_to_msg(expr1,expr2, msg) \
324 : do { \
325 : if (!((expr1) == (expr2))) { \
326 : libmesh_error_msg(std::setprecision(17) << "Assertion `" #expr1 " == " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl); \
327 : } } while (0)
328 :
329 : #define libmesh_assert_not_equal_to_msg(expr1,expr2, msg) \
330 : do { \
331 : if (!((expr1) != (expr2))) { \
332 : libmesh_error_msg(std::setprecision(17) << "Assertion `" #expr1 " != " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl); \
333 : } } while (0)
334 :
335 : template <template <class> class Comp>
336 : struct casting_compare {
337 :
338 : template <typename T1, typename T2>
339 963436770 : bool operator()(const T1 & e1, const T2 & e2) const
340 : {
341 : typedef typename std::decay<T1>::type DT1;
342 : typedef typename std::decay<T2>::type DT2;
343 1926873540 : return (Comp<DT2>()(static_cast<DT2>(e1), e2) &&
344 1926873540 : Comp<DT1>()(e1, static_cast<DT1>(e2)));
345 : }
346 :
347 : template <typename T1>
348 168626946 : bool operator()(const T1 & e1, const T1 & e2) const
349 : {
350 168626946 : return Comp<T1>()(e1, e2);
351 : }
352 : };
353 :
354 : #define libmesh_assert_less_msg(expr1,expr2, msg) \
355 : do { \
356 : if (!libMesh::casting_compare<std::less>()(expr1, expr2)) { \
357 : libmesh_error_msg(std::setprecision(17) << "Assertion `" #expr1 " < " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl); \
358 : } } while (0)
359 :
360 : #define libmesh_assert_greater_msg(expr1,expr2, msg) \
361 : do { \
362 : if (!libMesh::casting_compare<std::greater>()(expr1, expr2)) { \
363 : libmesh_error_msg(std::setprecision(17) << "Assertion `" #expr1 " > " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl); \
364 : } } while (0)
365 :
366 : #define libmesh_assert_less_equal_msg(expr1,expr2, msg) \
367 : do { \
368 : if (!libMesh::casting_compare<std::less_equal>()(expr1, expr2)) { \
369 : libmesh_error_msg(std::setprecision(17) << "Assertion `" #expr1 " <= " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl); \
370 : } } while (0)
371 :
372 : #define libmesh_assert_greater_equal_msg(expr1,expr2, msg) \
373 : do { \
374 : if (!libMesh::casting_compare<std::greater_equal>()(expr1, expr2)) { \
375 : libmesh_error_msg(std::setprecision(17) << "Assertion `" #expr1 " >= " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl); \
376 : } } while (0)
377 :
378 : #endif
379 :
380 :
381 : #define libmesh_assert(asserted) libmesh_assert_msg(asserted, "")
382 : #define libmesh_exceptionless_assert(asserted) libmesh_exceptionless_assert_msg(asserted, "")
383 : #define libmesh_assert_equal_to(expr1,expr2) libmesh_assert_equal_to_msg(expr1,expr2, "")
384 : #define libmesh_assert_not_equal_to(expr1,expr2) libmesh_assert_not_equal_to_msg(expr1,expr2, "")
385 : #define libmesh_assert_less(expr1,expr2) libmesh_assert_less_msg(expr1,expr2, "")
386 : #define libmesh_assert_greater(expr1,expr2) libmesh_assert_greater_msg(expr1,expr2, "")
387 : #define libmesh_assert_less_equal(expr1,expr2) libmesh_assert_less_equal_msg(expr1,expr2, "")
388 : #define libmesh_assert_greater_equal(expr1,expr2) libmesh_assert_greater_equal_msg(expr1,expr2, "")
389 :
390 : // The libmesh_error() macro prints a message and throws a LogicError
391 : // exception
392 : //
393 : // The libmesh_not_implemented() macro prints a message and throws a
394 : // NotImplemented exception
395 : //
396 : // The libmesh_file_error(const std::string & filename) macro prints a message
397 : // and throws a FileError exception
398 : //
399 : // The libmesh_convergence_failure() macro
400 : // throws a ConvergenceFailure exception
401 : #define libmesh_error_msg(msg) \
402 : do { \
403 : std::stringstream message_stream; \
404 : message_stream << msg << '\n'; \
405 : libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME, message_stream); \
406 : LIBMESH_THROW(libMesh::LogicError(message_stream.str())); \
407 : } while (0)
408 :
409 : #define libmesh_error() libmesh_error_msg("")
410 :
411 : #define libmesh_error_msg_if(cond, msg) \
412 : do { \
413 : if (cond) \
414 : libmesh_error_msg(msg); \
415 : } while (0)
416 :
417 : #define libmesh_exceptionless_error_msg(msg) \
418 : do { \
419 : libMesh::Threads::lock_singleton_spin_mutex(); \
420 : libMesh::err << msg << '\n'; \
421 : libMesh::Threads::unlock_singleton_spin_mutex(); \
422 : libmesh_try { libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); } \
423 : libmesh_catch (...) {} \
424 : std::terminate(); \
425 : } while (0)
426 :
427 : #define libmesh_exceptionless_error() libmesh_exceptionless_error_msg("")
428 :
429 : #define libmesh_not_implemented_msg(msg) \
430 : do { \
431 : std::stringstream message_stream; \
432 : message_stream << msg << '\n'; \
433 : libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME, message_stream); \
434 : LIBMESH_THROW(libMesh::NotImplemented(message_stream.str())); \
435 : } while (0)
436 :
437 : #define libmesh_not_implemented() libmesh_not_implemented_msg("")
438 :
439 : #define libmesh_file_error_msg(filename, msg) \
440 : do { \
441 : std::stringstream message_stream; \
442 : message_stream << msg << '\n'; \
443 : libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME, message_stream); \
444 : LIBMESH_THROW(libMesh::FileError(filename, message_stream.str())); \
445 : } while (0)
446 :
447 : #define libmesh_file_error(filename) libmesh_file_error_msg(filename,"")
448 :
449 : #define libmesh_convergence_failure() \
450 : do { \
451 : LIBMESH_THROW(libMesh::ConvergenceFailure()); \
452 : } while (0)
453 :
454 : // The libmesh_example_requires() macro prints a message and calls
455 : // "return 77;" if the condition specified by the macro is not true. This
456 : // macro is used in the example executables, which should run when the
457 : // configure-time libMesh options support them but which should exit
458 : // without failure otherwise.
459 : //
460 : // This macro only works in main(), because we have no better way than
461 : // "return" from main to immediately exit successfully - std::exit()
462 : // gets seen by at least some MPI stacks as failure.
463 : //
464 : // 77 is the automake code for a skipped test.
465 :
466 : #define libmesh_example_requires(condition, option) \
467 : do { \
468 : if (!(condition)) { \
469 : libMesh::out << "Configuring libMesh with " << option << " is required to run this example." << std::endl; \
470 : return 77; \
471 : } } while (0)
472 :
473 : // The libmesh_do_once macro helps us avoid redundant repeated
474 : // repetitions of the same warning messages
475 : #undef libmesh_do_once
476 : #define libmesh_do_once(do_this) \
477 : do { \
478 : static bool did_this_already = false; \
479 : if (!did_this_already) { \
480 : did_this_already = true; \
481 : do_this; \
482 : } } while (0)
483 :
484 :
485 : // The libmesh_warning macro outputs a file/line/time stamped warning
486 : // message, if warnings are enabled.
487 : #ifdef LIBMESH_ENABLE_WARNINGS
488 : #define libmesh_warning(message) \
489 : libmesh_do_once(libMesh::out << message \
490 : << __FILE__ << ", line " << __LINE__ << ", compiled " << LIBMESH_DATE << " at " << LIBMESH_TIME << " ***" << std::endl;)
491 : #else
492 : #define libmesh_warning(message) ((void) 0)
493 : #endif
494 :
495 : // The libmesh_experimental macro warns that you are using
496 : // bleeding-edge code
497 : #undef libmesh_experimental
498 : #define libmesh_experimental() \
499 : libmesh_warning("*** Warning, This code is untested, experimental, or likely to see future API changes: ");
500 :
501 :
502 : // The libmesh_deprecated macro warns that you are using obsoleted code
503 : #undef libmesh_deprecated
504 : #ifndef LIBMESH_ENABLE_DEPRECATED
505 : #define libmesh_deprecated() \
506 : libmesh_error_msg("*** Error, This code is deprecated, and likely to be removed in future library versions! ");
507 : #else
508 : #define libmesh_deprecated() \
509 : libmesh_warning("*** Warning, This code is deprecated, and likely to be removed in future library versions! ");
510 : #endif
511 :
512 : // A function template for ignoring unused variables. This is a way
513 : // to shut up unused variable compiler warnings on a case by case
514 : // basis.
515 303935873 : template<class ...Args> inline void libmesh_ignore( const Args&... ) { }
516 :
517 :
518 : // A workaround for the lack of C++17 merge() support in some
519 : // compilers
520 :
521 : #ifdef LIBMESH_HAVE_CXX17_SPLICING
522 : template <typename T>
523 47976 : void libmesh_merge_move(T & target, T & source)
524 : {
525 47976 : target.merge(std::move(source));
526 1420728 : }
527 : #else
528 : template <typename T>
529 : void libmesh_merge_move(T & target, T & source)
530 : {
531 : target.insert(source.begin(), source.end());
532 : source.clear(); // Avoid forwards-incompatibility
533 : }
534 : #endif // LIBMESH_HAVE_CXX17_SPLICING
535 :
536 : /**
537 : * Mostly system independent demangler
538 : */
539 : std::string demangle(const char * name);
540 :
541 : // cast_ref and cast_ptr do a dynamic cast and assert
542 : // the result, if we have RTTI enabled and we're in debug or
543 : // development modes, but they just do a faster static cast if we're
544 : // in optimized mode.
545 : //
546 : // Use these casts when you're certain that a cast will succeed in
547 : // correct code but you want to be able to double-check.
548 : template <typename Tnew, typename Told>
549 8584090 : inline Tnew cast_ref(Told & oldvar)
550 : {
551 : #if !defined(NDEBUG) && defined(LIBMESH_HAVE_RTTI) && defined(LIBMESH_ENABLE_EXCEPTIONS)
552 : try
553 : {
554 8584090 : Tnew newvar = dynamic_cast<Tnew>(oldvar);
555 8584090 : return newvar;
556 : }
557 0 : catch (std::bad_cast &)
558 : {
559 0 : libMesh::err << "Failed to convert " << demangle(typeid(Told).name())
560 0 : << " reference to " << demangle(typeid(Tnew).name())
561 0 : << std::endl;
562 0 : libMesh::err << "The " << demangle(typeid(Told).name())
563 0 : << " appears to be a "
564 0 : << demangle(typeid(*(&oldvar)).name()) << std::endl;
565 0 : libmesh_error();
566 : }
567 : #else
568 : return(static_cast<Tnew>(oldvar));
569 : #endif
570 : }
571 :
572 : // We use two different function names to avoid an odd overloading
573 : // ambiguity bug with icc 10.1.008
574 : template <typename Tnew, typename Told>
575 745469264 : inline Tnew cast_ptr (Told * oldvar)
576 : {
577 : #if !defined(NDEBUG) && defined(LIBMESH_HAVE_RTTI)
578 745469264 : Tnew newvar = dynamic_cast<Tnew>(oldvar);
579 745469264 : if (!newvar)
580 : {
581 0 : libMesh::err << "Failed to convert " << demangle(typeid(Told).name())
582 0 : << " pointer to " << demangle(typeid(Tnew).name())
583 0 : << std::endl;
584 0 : libMesh::err << "The " << demangle(typeid(Told).name())
585 0 : << " appears to be a "
586 0 : << demangle(typeid(*oldvar).name()) << std::endl;
587 0 : libmesh_error();
588 : }
589 745469264 : return newvar;
590 : #else
591 : return(static_cast<Tnew>(oldvar));
592 : #endif
593 : }
594 :
595 :
596 : #ifdef LIBMESH_ENABLE_DEPRECATED
597 : template <typename Tnew, typename Told>
598 : inline Tnew libmesh_cast_ptr (Told * oldvar)
599 : {
600 : libmesh_deprecated();
601 :
602 : // we use the less redundantly named libMesh::cast_ptr now
603 : return cast_ptr<Tnew>(oldvar);
604 : }
605 : #endif // LIBMESH_ENABLE_DEPRECATED
606 :
607 :
608 : // cast_int asserts that the value of the castee is within the
609 : // bounds which are exactly representable by the output type, if we're
610 : // in debug or development modes, but it just does a faster static
611 : // cast if we're in optimized mode.
612 : //
613 : // Use these casts when you're certain that a cast will succeed in
614 : // correct code but you want to be able to double-check.
615 : template <typename Tnew, typename Told>
616 8261421010 : inline Tnew cast_int (Told oldvar)
617 : {
618 133601451 : libmesh_assert_equal_to
619 : (oldvar, static_cast<Told>(static_cast<Tnew>(oldvar)));
620 :
621 20851184189 : return(static_cast<Tnew>(oldvar));
622 : }
623 :
624 :
625 : template <typename Tnew, typename Told>
626 : inline Tnew libmesh_cast_int (Told oldvar)
627 : {
628 : // we use the less redundantly named libMesh::cast_int now
629 : return cast_int<Tnew>(oldvar);
630 : }
631 :
632 : /**
633 : * This is a helper variable template for cases when we want to use a default compile-time
634 : * error with constexpr-based if conditions. The templating delays the triggering
635 : * of the static assertion until the template is instantiated.
636 : */
637 : template <class T>
638 : constexpr std::false_type always_false{};
639 :
640 : static constexpr std::size_t libmesh_dim = LIBMESH_DIM;
641 :
642 : // build a integer representation of version
643 : #define LIBMESH_VERSION_ID(major,minor,patch) (((major) << 16) | ((minor) << 8) | ((patch) & 0xFF))
644 :
645 :
646 : // libmesh_override is simply a synonym for override as we now require
647 : // a C++11 compiler that supports this keyword.
648 : #define libmesh_override override
649 :
650 : // libmesh_delete is simply a synonym for '=delete' as we now require
651 : // a C++11 compiler that supports this keyword.
652 : #define libmesh_delete =delete
653 :
654 : // libmesh_final is simply a synonym for 'final' as we now require
655 : // a C++11 compiler that supports this keyword.
656 : #define libmesh_final final
657 :
658 : // Define backwards-compatible fallthrough attribute. We could
659 : // eventually also add support for other compiler-specific fallthrough
660 : // attributes.
661 : #ifdef LIBMESH_HAVE_CXX17_FALLTHROUGH_ATTRIBUTE
662 : #define libmesh_fallthrough() [[fallthrough]]
663 : #elif defined(LIBMESH_HAVE_DOUBLE_UNDERSCORE_ATTRIBUTE_FALLTHROUGH)
664 : #define libmesh_fallthrough() __attribute__((fallthrough))
665 : #else
666 : #define libmesh_fallthrough() ((void) 0)
667 : #endif
668 :
669 : } // namespace libMesh
670 :
671 :
672 : // Backwards compatibility
673 : namespace libMeshEnums
674 : {
675 : using namespace libMesh;
676 : }
677 :
678 : // Backwards compatibility with pre-TIMPI reference
679 : namespace TIMPI {}
680 :
681 : namespace libMesh {
682 : namespace Parallel {
683 : using namespace TIMPI;
684 : }
685 : }
686 :
687 :
688 : // Here we add missing types to the standard namespace. For example,
689 : // std::max(double, float) etc... are well behaved but not defined
690 : // by the standard. This also includes workarounds for super-strict
691 : // implementations, for example Sun Studio and PGI C++. However,
692 : // this necessarily requires breaking the ISO-C++ standard, and is
693 : // really just a hack. As such, only do it if we are building the
694 : // libmesh library itself. Specifically, *DO NOT* export this to
695 : // user code or install this header.
696 : //
697 : // We put this at the end of libmesh_common.h so we can make use of
698 : // any exotic definitions of Real above.
699 : #ifdef LIBMESH_IS_COMPILING_ITSELF
700 : # include "libmesh/libmesh_augment_std_namespace.h"
701 : #endif
702 :
703 :
704 : #ifdef _MSC_VER
705 : #define LIBMESH_EXPORT __declspec(dllexport)
706 : #else
707 : #define LIBMESH_EXPORT
708 : #endif
709 :
710 :
711 : #endif // LIBMESH_LIBMESH_COMMON_H
|