libMesh
libmesh_common.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 // The library configuration options
24 #include "libmesh/libmesh_config.h"
25 
26 // Use actual timestamps or constant dummies (to aid ccache)
27 #ifdef LIBMESH_ENABLE_TIMESTAMPS
28 # define LIBMESH_TIME __TIME__
29 # define LIBMESH_DATE __DATE__
30 #else
31 # define LIBMESH_TIME "notime"
32 # define LIBMESH_DATE "nodate"
33 #endif
34 
35 // C/C++ includes everyone should know about
36 #include <cstdlib>
37 #ifdef __PGI
38 // BSK, Thu Feb 20 08:32:06 CST 2014 - For some reason, unless PGI gets
39 // <cmath> early this nonsense shows up:
40 // "/software/x86_64/pgi/12.9/linux86-64/12.9/include/CC/cmath", line 57: error:
41 // the global scope has no "abs"
42 // using _STLP_VENDOR_CSTD::abs;
43 // So include <cmath> as early as possible under the PGI compilers.
44 # include <cmath>
45 #endif
46 #include <complex>
47 #include <typeinfo> // std::bad_cast
48 
49 
50 // Include the MPI definition
51 #ifdef LIBMESH_HAVE_MPI
52 # include "libmesh/ignore_warnings.h"
53 # include <mpi.h>
54 # include "libmesh/restore_warnings.h"
55 #endif
56 
57 // _basic_ library functionality
58 #include "libmesh/libmesh_base.h"
59 #include "libmesh/libmesh_exceptions.h"
60 
61 // Proxy class for libMesh::out/err output
62 #include "libmesh/ostream_proxy.h"
63 
64 // Here we add missing types to the standard namespace. For example,
65 // std::max(double, float) etc... are well behaved but not defined
66 // by the standard. This also includes workarounds for super-strict
67 // implementations, for example Sun Studio and PGI C++. However,
68 // this necessarily requires breaking the ISO-C++ standard, and is
69 // really just a hack. As such, only do it if we are building the
70 // libmesh library itself. Specifically, *DO NOT* export this to
71 // user code or install this header.
72 #ifdef LIBMESH_IS_COMPILING_ITSELF
73 # include "libmesh/libmesh_augment_std_namespace.h"
74 #endif
75 
76 // Make sure the libmesh_nullptr define is available for backwards
77 // compatibility, although we no longer use it in the library.
78 #include "libmesh/libmesh_nullptr.h"
79 
80 namespace libMesh
81 {
82 
83 
84 // A namespace for functions used in the bodies of the macros below.
85 // The macros generally call these functions with __FILE__, __LINE__,
86 // __DATE__, and __TIME__ in the appropriate order. These should not
87 // be called by users directly! The implementations can be found in
88 // libmesh_common.C.
89 namespace MacroFunctions
90 {
91 void here(const char * file, int line, const char * date, const char * time);
92 void stop(const char * file, int line, const char * date, const char * time);
93 void report_error(const char * file, int line, const char * date, const char * time);
94 }
95 
96 // Undefine any existing macros
97 #ifdef Real
98 # undef Real
99 #endif
100 
101 //#ifdef REAL
102 //# undef REAL
103 //#endif
104 
105 #ifdef Complex
106 # undef Complex
107 #endif
108 
109 #ifdef COMPLEX
110 # undef COMPLEX
111 #endif
112 
113 #ifdef MPI_REAL
114 # undef MPI_REAL
115 #endif
116 
117 // Check to see if TOLERANCE has been defined by another
118 // package, if so we might want to change the name...
119 #ifdef TOLERANCE
120 DIE A HORRIBLE DEATH HERE...
121 # undef TOLERANCE
122 #endif
123 
124 
125 
126 // Define the type to use for real numbers
127 
128 typedef LIBMESH_DEFAULT_SCALAR_TYPE Real;
129 
130 // Define a corresponding tolerance. This is what should be
131 // considered "good enough" when doing floating point comparisons.
132 // For example, v == 0 is changed to std::abs(v) < TOLERANCE.
133 
134 #ifndef LIBMESH_DEFAULT_SINGLE_PRECISION
135 #ifdef LIBMESH_DEFAULT_TRIPLE_PRECISION
136 static const Real TOLERANCE = 1.e-8;
137 # define MPI_REAL MPI_LONG_DOUBLE
138 #else
139 static const Real TOLERANCE = 1.e-6;
140 # define MPI_REAL MPI_DOUBLE
141 #endif
142 #else
143 static const Real TOLERANCE = 2.5e-3;
144 # define MPI_REAL MPI_FLOAT
145 #endif
146 
147 // Define the type to use for complex numbers
148 // Always use std::complex<double>, as required by Petsc?
149 // If your version of Petsc doesn't support
150 // std::complex<other_precision>, then you'd better just leave
151 // Real==double
152 typedef std::complex<Real> Complex;
153 typedef std::complex<Real> COMPLEX;
154 
155 
156 // Helper functions for complex/real numbers
157 // to clean up #ifdef LIBMESH_USE_COMPLEX_NUMBERS elsewhere
158 template<typename T> inline T libmesh_real(T a) { return a; }
159 template<typename T> inline T libmesh_conj(T a) { return a; }
160 
161 template<typename T>
162 inline T libmesh_real(std::complex<T> a) { return std::real(a); }
163 
164 template<typename T>
165 inline std::complex<T> libmesh_conj(std::complex<T> a) { return std::conj(a); }
166 
167 // std::isnan() is in <cmath> as of C++11.
168 template <typename T>
169 inline bool libmesh_isnan(T x) { return std::isnan(x); }
170 
171 template <typename T>
172 inline bool libmesh_isnan(std::complex<T> a)
173 { return (std::isnan(std::real(a)) || std::isnan(std::imag(a))); }
174 
175 // std::isinf() is in <cmath> as of C++11.
176 template <typename T>
177 inline bool libmesh_isinf(T x) { return std::isinf(x); }
178 
179 template <typename T>
180 inline bool libmesh_isinf(std::complex<T> a)
181 { return (std::isinf(std::real(a)) || std::isinf(std::imag(a))); }
182 
183 // Define the value type for unknowns in simulations.
184 // This is either Real or Complex, depending on how
185 // the library was configures
186 #if defined (LIBMESH_USE_REAL_NUMBERS)
187 typedef Real Number;
188 #elif defined (LIBMESH_USE_COMPLEX_NUMBERS)
189 typedef Complex Number;
190 #else
191 DIE A HORRIBLE DEATH HERE...
192 #endif
193 
194 
195 // Define the value type for error estimates.
196 // Since AMR/C decisions don't have to be precise,
197 // we default to float for memory efficiency.
198 typedef float ErrorVectorReal;
199 #define MPI_ERRORVECTORREAL MPI_FLOAT
200 
201 
202 #ifdef LIBMESH_HAVE_MPI
203 
207 extern MPI_Comm GLOBAL_COMM_WORLD;
208 #else
209 
214 extern int GLOBAL_COMM_WORLD;
215 #endif
216 
217 // Let's define a couple output streams - these will default
218 // to cout/cerr, but LibMeshInit (or the user) can also set them to
219 // something more sophisticated.
220 //
221 // We use a proxy class rather than references so they can be
222 // reseated at runtime.
223 
224 extern OStreamProxy out;
225 extern OStreamProxy err;
226 
227 // This global variable is to help us deprecate AutoPtr. We can't
228 // just use libmesh_deprecated() because then you get one print out
229 // per template instantiation, instead of one total print out.
230 extern bool warned_about_auto_ptr;
231 
232 // These are useful macros that behave like functions in the code.
233 // If you want to make sure you are accessing a section of code just
234 // stick a libmesh_here(); in it, for example
235 #define libmesh_here() \
236  do { \
237  libMesh::MacroFunctions::here(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
238  } while (0)
239 
240 // the libmesh_stop() macro will stop the code until a SIGCONT signal
241 // is received. This is useful, for example, when determining the
242 // memory used by a given operation. A libmesh_stop() could be
243 // inserted before and after a questionable operation and the delta
244 // memory can be obtained from a ps or top. This macro only works for
245 // serial cases.
246 #define libmesh_stop() \
247  do { \
248  libMesh::MacroFunctions::stop(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
249  } while (0)
250 
251 // The libmesh_dbg_var() macro indicates that an argument to a function
252 // is used only in debug mode (i.e., when NDEBUG is not defined).
253 #ifndef NDEBUG
254 #define libmesh_dbg_var(var) var
255 #else
256 #define libmesh_dbg_var(var)
257 #endif
258 
259 // The libmesh_assert() macro acts like C's assert(), but throws a
260 // libmesh_error() (including stack trace, etc) instead of just exiting
261 #ifdef NDEBUG
262 
263 #define libmesh_assert_msg(asserted, msg) ((void) 0)
264 #define libmesh_exceptionless_assert_msg(asserted, msg) ((void) 0)
265 #define libmesh_assert_equal_to_msg(expr1,expr2, msg) ((void) 0)
266 #define libmesh_assert_not_equal_to_msg(expr1,expr2, msg) ((void) 0)
267 #define libmesh_assert_less_msg(expr1,expr2, msg) ((void) 0)
268 #define libmesh_assert_greater_msg(expr1,expr2, msg) ((void) 0)
269 #define libmesh_assert_less_equal_msg(expr1,expr2, msg) ((void) 0)
270 #define libmesh_assert_greater_equal_msg(expr1,expr2, msg) ((void) 0)
271 
272 #else
273 
274 #define libmesh_assert_msg(asserted, msg) \
275  do { \
276  if (!(asserted)) { \
277  libMesh::err << "Assertion `" #asserted "' failed." << std::endl; \
278  libmesh_error_msg(msg); \
279  } } while (0)
280 
281 #define libmesh_exceptionless_assert_msg(asserted, msg) \
282  do { \
283  if (!(asserted)) { \
284  libMesh::err << "Assertion `" #asserted "' failed." << std::endl; \
285  libmesh_exceptionless_error(); \
286  } } while (0)
287 
288 #define libmesh_assert_equal_to_msg(expr1,expr2, msg) \
289  do { \
290  if (!(expr1 == expr2)) { \
291  libMesh::err << "Assertion `" #expr1 " == " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl; \
292  libmesh_error(); \
293  } } while (0)
294 
295 #define libmesh_assert_not_equal_to_msg(expr1,expr2, msg) \
296  do { \
297  if (!(expr1 != expr2)) { \
298  libMesh::err << "Assertion `" #expr1 " != " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl; \
299  libmesh_error(); \
300  } } while (0)
301 
302 #define libmesh_assert_less_msg(expr1,expr2, msg) \
303  do { \
304  if (!(expr1 < expr2)) { \
305  libMesh::err << "Assertion `" #expr1 " < " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl; \
306  libmesh_error(); \
307  } } while (0)
308 
309 #define libmesh_assert_greater_msg(expr1,expr2, msg) \
310  do { \
311  if (!(expr1 > expr2)) { \
312  libMesh::err << "Assertion `" #expr1 " > " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl; \
313  libmesh_error(); \
314  } } while (0)
315 
316 #define libmesh_assert_less_equal_msg(expr1,expr2, msg) \
317  do { \
318  if (!(expr1 <= expr2)) { \
319  libMesh::err << "Assertion `" #expr1 " <= " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl; \
320  libmesh_error(); \
321  } } while (0)
322 
323 #define libmesh_assert_greater_equal_msg(expr1,expr2, msg) \
324  do { \
325  if (!(expr1 >= expr2)) { \
326  libMesh::err << "Assertion `" #expr1 " >= " #expr2 "' failed.\n" #expr1 " = " << (expr1) << "\n" #expr2 " = " << (expr2) << '\n' << msg << std::endl; \
327  libmesh_error(); \
328  } } while (0)
329 #endif
330 
331 
332 #define libmesh_assert(asserted) libmesh_assert_msg(asserted, "")
333 #define libmesh_exceptionless_assert(asserted) libmesh_exceptionless_assert_msg(asserted, "")
334 #define libmesh_assert_equal_to(expr1,expr2) libmesh_assert_equal_to_msg(expr1,expr2, "")
335 #define libmesh_assert_not_equal_to(expr1,expr2) libmesh_assert_not_equal_to_msg(expr1,expr2, "")
336 #define libmesh_assert_less(expr1,expr2) libmesh_assert_less_msg(expr1,expr2, "")
337 #define libmesh_assert_greater(expr1,expr2) libmesh_assert_greater_msg(expr1,expr2, "")
338 #define libmesh_assert_less_equal(expr1,expr2) libmesh_assert_less_equal_msg(expr1,expr2, "")
339 #define libmesh_assert_greater_equal(expr1,expr2) libmesh_assert_greater_equal_msg(expr1,expr2, "")
340 
341 // The libmesh_error() macro prints a message and throws a LogicError
342 // exception
343 //
344 // The libmesh_not_implemented() macro prints a message and throws a
345 // NotImplemented exception
346 //
347 // The libmesh_file_error(const std::string & filename) macro prints a message
348 // and throws a FileError exception
349 //
350 // The libmesh_convergence_failure() macro
351 // throws a ConvergenceFailure exception
352 #define libmesh_error_msg(msg) \
353  do { \
354  libMesh::err << msg << std::endl; \
355  std::stringstream msg_stream; \
356  msg_stream << msg; \
357  libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
358  LIBMESH_THROW(libMesh::LogicError(msg_stream.str())); \
359  } while (0)
360 
361 #define libmesh_error() libmesh_error_msg("")
362 
363 #define libmesh_exceptionless_error_msg(msg) \
364  do { \
365  libMesh::err << msg << std::endl; \
366  libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
367  std::terminate(); \
368  } while (0)
369 
370 #define libmesh_exceptionless_error() libmesh_exceptionless_error_msg("")
371 
372 #define libmesh_not_implemented_msg(msg) \
373  do { \
374  libMesh::err << msg << std::endl; \
375  libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
376  LIBMESH_THROW(libMesh::NotImplemented()); \
377  } while (0)
378 
379 #define libmesh_not_implemented() libmesh_not_implemented_msg("")
380 
381 #define libmesh_file_error_msg(filename, msg) \
382  do { \
383  libMesh::err << "Error with file `" << filename << "'" << std::endl; \
384  libMesh::MacroFunctions::report_error(__FILE__, __LINE__, LIBMESH_DATE, LIBMESH_TIME); \
385  libMesh::err << msg << std::endl; \
386  LIBMESH_THROW(libMesh::FileError(filename)); \
387  } while (0)
388 
389 #define libmesh_file_error(filename) libmesh_file_error_msg(filename,"")
390 
391 #define libmesh_convergence_failure() \
392  do { \
393  LIBMESH_THROW(libMesh::ConvergenceFailure()); \
394  } while (0)
395 
396 // The libmesh_example_requires() macro prints a message and calls
397 // "return 77;" if the condition specified by the macro is not true. This
398 // macro is used in the example executables, which should run when the
399 // configure-time libMesh options support them but which should exit
400 // without failure otherwise.
401 //
402 // This macro only works in main(), because we have no better way than
403 // "return" from main to immediately exit successfully - std::exit()
404 // gets seen by at least some MPI stacks as failure.
405 //
406 // 77 is the automake code for a skipped test.
407 
408 #define libmesh_example_requires(condition, option) \
409  do { \
410  if (!(condition)) { \
411  libMesh::out << "Configuring libMesh with " << option << " is required to run this example." << std::endl; \
412  return 77; \
413  } } while (0)
414 
415 // The libmesh_do_once macro helps us avoid redundant repeated
416 // repetitions of the same warning messages
417 #undef libmesh_do_once
418 #define libmesh_do_once(do_this) \
419  do { \
420  static bool did_this_already = false; \
421  if (!did_this_already) { \
422  did_this_already = true; \
423  do_this; \
424  } } while (0)
425 
426 
427 // The libmesh_warning macro outputs a file/line/time stamped warning
428 // message, if warnings are enabled.
429 #ifdef LIBMESH_ENABLE_WARNINGS
430 #define libmesh_warning(message) \
431  libmesh_do_once(libMesh::out << message \
432  << __FILE__ << ", line " << __LINE__ << ", compiled " << LIBMESH_DATE << " at " << LIBMESH_TIME << " ***" << std::endl;)
433 #else
434 #define libmesh_warning(message) ((void) 0)
435 #endif
436 
437 // The libmesh_experimental macro warns that you are using
438 // bleeding-edge code
439 #undef libmesh_experimental
440 #define libmesh_experimental() \
441  libmesh_warning("*** Warning, This code is untested, experimental, or likely to see future API changes: ");
442 
443 
444 // The libmesh_deprecated macro warns that you are using obsoleted code
445 #undef libmesh_deprecated
446 #ifndef LIBMESH_ENABLE_DEPRECATED
447 #define libmesh_deprecated() \
448  libmesh_error_msg("*** Error, This code is deprecated, and likely to be removed in future library versions! ");
449 #else
450 #define libmesh_deprecated() \
451  libmesh_warning("*** Warning, This code is deprecated, and likely to be removed in future library versions! ");
452 #endif
453 
454 // A function template for ignoring unused variables. This is a way
455 // to shut up unused variable compiler warnings on a case by case
456 // basis.
457 template<class ...Args> inline void libmesh_ignore( const Args&... ) { }
458 
459 
460 // cast_ref and cast_ptr do a dynamic cast and assert
461 // the result, if we have RTTI enabled and we're in debug or
462 // development modes, but they just do a faster static cast if we're
463 // in optimized mode.
464 //
465 // Use these casts when you're certain that a cast will succeed in
466 // correct code but you want to be able to double-check.
467 template <typename Tnew, typename Told>
468 inline Tnew cast_ref(Told & oldvar)
469 {
470 #if !defined(NDEBUG) && defined(LIBMESH_HAVE_RTTI) && defined(LIBMESH_ENABLE_EXCEPTIONS)
471  try
472  {
473  Tnew newvar = dynamic_cast<Tnew>(oldvar);
474  return newvar;
475  }
476  catch (std::bad_cast &)
477  {
478  libMesh::err << "Failed to convert " << typeid(Told).name()
479  << " reference to " << typeid(Tnew).name()
480  << std::endl;
481  libMesh::err << "The " << typeid(Told).name()
482  << " appears to be a "
483  << typeid(*(&oldvar)).name() << std::endl;
484  libmesh_error();
485  }
486 #else
487  return(static_cast<Tnew>(oldvar));
488 #endif
489 }
490 
491 #ifdef LIBMESH_ENABLE_DEPRECATED
492 template <typename Tnew, typename Told>
493 inline Tnew libmesh_cast_ref(Told & oldvar)
494 {
495  // we use the less redundantly named libMesh::cast_ref now
496  libmesh_deprecated();
497  return cast_ref<Tnew>(oldvar);
498 }
499 #endif
500 
501 // We use two different function names to avoid an odd overloading
502 // ambiguity bug with icc 10.1.008
503 template <typename Tnew, typename Told>
504 inline Tnew cast_ptr (Told * oldvar)
505 {
506 #if !defined(NDEBUG) && defined(LIBMESH_HAVE_RTTI)
507  Tnew newvar = dynamic_cast<Tnew>(oldvar);
508  if (!newvar)
509  {
510  libMesh::err << "Failed to convert " << typeid(Told).name()
511  << " pointer to " << typeid(Tnew).name()
512  << std::endl;
513  libMesh::err << "The " << typeid(Told).name()
514  << " appears to be a "
515  << typeid(*oldvar).name() << std::endl;
516  libmesh_error();
517  }
518  return newvar;
519 #else
520  return(static_cast<Tnew>(oldvar));
521 #endif
522 }
523 
524 
525 template <typename Tnew, typename Told>
526 inline Tnew libmesh_cast_ptr (Told * oldvar)
527 {
528  // we use the less redundantly named libMesh::cast_ptr now
529  return cast_ptr<Tnew>(oldvar);
530 }
531 
532 
533 // cast_int asserts that the value of the castee is within the
534 // bounds which are exactly representable by the output type, if we're
535 // in debug or development modes, but it just does a faster static
536 // cast if we're in optimized mode.
537 //
538 // Use these casts when you're certain that a cast will succeed in
539 // correct code but you want to be able to double-check.
540 template <typename Tnew, typename Told>
541 inline Tnew cast_int (Told oldvar)
542 {
543  libmesh_assert_equal_to
544  (oldvar, static_cast<Told>(static_cast<Tnew>(oldvar)));
545 
546  return(static_cast<Tnew>(oldvar));
547 }
548 
549 
550 template <typename Tnew, typename Told>
551 inline Tnew libmesh_cast_int (Told oldvar)
552 {
553  // we use the less redundantly named libMesh::cast_int now
554  return cast_int<Tnew>(oldvar);
555 }
556 
557 
558 // build a integer representation of version
559 #define LIBMESH_VERSION_ID(major,minor,patch) (((major) << 16) | ((minor) << 8) | ((patch) & 0xFF))
560 
561 
562 // libmesh_override is simply a synonym for override as we now require
563 // a C++11 compiler that supports this keyword.
564 #define libmesh_override override
565 
566 // libmesh_delete is simply a synonym for '=delete' as we now require
567 // a C++11 compiler that supports this keyword.
568 #define libmesh_delete =delete
569 
570 // libmesh_final is simply a synonym for 'final' as we now require
571 // a C++11 compiler that supports this keyword.
572 #define libmesh_final final
573 
574 // Define backwards-compatible fallthrough attribute. We could
575 // eventually also add support for other compiler-specific fallthrough
576 // attributes.
577 #ifdef LIBMESH_HAVE_CXX17_FALLTHROUGH_ATTRIBUTE
578 #define libmesh_fallthrough() [[fallthrough]]
579 #elif defined(LIBMESH_HAVE_DOUBLE_UNDERSCORE_ATTRIBUTE_FALLTHROUGH)
580 #define libmesh_fallthrough() __attribute__((fallthrough))
581 #else
582 #define libmesh_fallthrough() ((void) 0)
583 #endif
584 
585 } // namespace libMesh
586 
587 
588 // Backwards compatibility
589 namespace libMeshEnums
590 {
591 using namespace libMesh;
592 }
593 
594 #endif // LIBMESH_LIBMESH_COMMON_H
T libmesh_real(T a)
std::string name(const ElemQuality q)
T libmesh_conj(T a)
Tnew cast_ref(Told &oldvar)
Tnew cast_ptr(Told *oldvar)
bool warned_about_auto_ptr
static const Real TOLERANCE
MPI_Comm GLOBAL_COMM_WORLD
MPI Communicator used to initialize libMesh.
Tnew libmesh_cast_ptr(Told *oldvar)
The libMesh namespace provides an interface to certain functionality in the library.
bool libmesh_isnan(T x)
DIE A HORRIBLE DEATH HERE typedef float ErrorVectorReal
bool libmesh_isinf(T x)
Tnew cast_int(Told oldvar)
This class is intended to be reseatable like a pointer-to-ostream for flexibility, but to look like a reference when used to produce less awkward user code.
Definition: ostream_proxy.h:42
void libmesh_ignore(const Args &...)
std::complex< Real > COMPLEX
void stop(const char *file, int line, const char *date, const char *time)
Tnew libmesh_cast_ref(Told &oldvar)
OStreamProxy out
std::complex< Real > Complex
OStreamProxy err
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void here(const char *file, int line, const char *date, const char *time)
void report_error(const char *file, int line, const char *date, const char *time)
Tnew libmesh_cast_int(Told oldvar)