Reference documentation for deal.II version 9.2.0
|
This module contains classes that are used in the exception mechanism of deal.II. More...
Modules | |
TrilinosWrappers | |
The classes in this module are wrappers around functionality provided by the Trilinos library. They provide a modern object-oriented interface that is compatible with the interfaces of the other linear algebra classes in deal.II. All classes and functions in this group reside in a namespace TrilinosWrappers . | |
Macros | |
#define | DeclException0(Exception0) |
#define | DeclExceptionMsg(Exception, defaulttext) |
#define | DeclException1(Exception1, type1, outsequence) |
#define | DeclException2(Exception2, type1, type2, outsequence) |
#define | DeclException3(Exception3, type1, type2, type3, outsequence) |
#define | DeclException4(Exception4, type1, type2, type3, type4, outsequence) |
#define | DeclException5( Exception5, type1, type2, type3, type4, type5, outsequence) |
#define | Assert(cond, exc) |
#define | AssertNothrow(cond, exc) |
#define | AssertThrow(cond, exc) |
#define | AssertDimension(dim1, dim2) |
#define | AssertVectorVectorDimension(VEC, DIM1, DIM2) |
#define | AssertIndexRange(index, range) |
#define | AssertIsFinite(number) |
#define | AssertThrowMPI(error_code) AssertThrow(error_code == MPI_SUCCESS, ::ExcMPI(error_code)) |
#define | AssertCuda(error_code) |
#define | AssertNothrowCuda(error_code) |
#define | AssertCudaKernel() |
#define | AssertCusparse(error_code) |
#define | AssertNothrowCusparse(error_code) |
#define | AssertCusolver(error_code) |
Data storage and access | |
static ::ExceptionBase & | GeneralDataStorage::ExcNameNotFound (std::string arg1) |
static ::ExceptionBase & | GeneralDataStorage::ExcNameHasBeenFound (std::string arg1) |
static ::ExceptionBase & | GeneralDataStorage::ExcTypeMismatch (std::string arg1, const char *arg2, const char *arg3) |
Bitwise operators | |
std::ostream & | Differentiation::SD::operator<< (std::ostream &stream, const Expression &expression) |
std::istream & | Differentiation::SD::operator>> (std::istream &stream, Expression &expression) |
Comparison operators | |
Expression | Differentiation::SD::operator== (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator!= (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator< (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator> (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator<= (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator>= (const Expression &lhs, const Expression &rhs) |
Logical operators | |
Expression | Differentiation::SD::operator! (const Expression &expression) |
Expression | Differentiation::SD::operator& (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator| (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator^ (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator&& (const Expression &lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator|| (const Expression &lhs, const Expression &rhs) |
Mathematical operators | |
Expression | Differentiation::SD::operator+ (Expression lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator- (Expression lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator* (Expression lhs, const Expression &rhs) |
Expression | Differentiation::SD::operator/ (Expression lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator+ (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator+ (const Expression &lhs, const NumberType &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator- (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator- (const Expression &lhs, const NumberType &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator* (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator* (const Expression &lhs, const NumberType &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator/ (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | Differentiation::SD::operator/ (const Expression &lhs, const NumberType &rhs) |
Constructors | |
Differentiation::SD::Expression::Expression () | |
Differentiation::SD::Expression::Expression (const bool value) | |
template<typename NumberType , typename = typename std::enable_if< std::is_arithmetic<NumberType>::value>::type> | |
Differentiation::SD::Expression::Expression (const NumberType &value) | |
template<typename NumberType , typename = typename std::enable_if< std::is_arithmetic<NumberType>::value>::type> | |
Differentiation::SD::Expression::Expression (const std::complex< NumberType > &value) | |
Differentiation::SD::Expression::Expression (const SymEngine::integer_class &value) | |
template<typename NumberType , typename = typename std::enable_if< std::is_integral<NumberType>::value>::type> | |
Differentiation::SD::Expression::Expression (const NumberType &numerator, const NumberType &denominator) | |
Differentiation::SD::Expression::Expression (const SymEngine::rational_class &value) | |
Differentiation::SD::Expression::Expression (const Expression &condition, const Expression &expression_if_true, const Expression &expression_if_false) | |
Differentiation::SD::Expression::Expression (const std::vector< std::pair< Expression, Expression >> &condition_expression, const Expression &expression_otherwise) | |
Differentiation::SD::Expression::Expression (const std::vector< std::pair< Expression, Expression >> &condition_expression) | |
Differentiation::SD::Expression::Expression (const char *symbol) | |
Differentiation::SD::Expression::Expression (const std::string &symb_expr, const bool parse_as_expression=false) | |
Differentiation::SD::Expression::Expression (const std::string &symbol_func, const types::symbol_vector &arguments) | |
Differentiation::SD::Expression::Expression (const Expression &rhs)=default | |
Differentiation::SD::Expression::Expression (const SymEngine::Expression &rhs) | |
Differentiation::SD::Expression::Expression (const SymEngine::RCP< const SymEngine::Basic > &rhs) | |
Differentiation::SD::Expression::Expression (Expression &&rhs)=default | |
Differentiation::SD::Expression::Expression (SymEngine::RCP< const SymEngine::Basic > &&rhs) | |
virtual | Differentiation::SD::Expression::~Expression ()=default |
Expression & | Differentiation::SD::Expression::parse (const std::string &expression) |
std::ostream & | Differentiation::SD::Expression::print (std::ostream &stream) const |
void | Differentiation::SD::Expression::save (std::ostream &stream) const |
void | Differentiation::SD::Expression::load (std::istream &stream) |
template<class Archive > | |
void | Differentiation::SD::Expression::save (Archive &archive, const unsigned int version) const |
template<class Archive > | |
void | Differentiation::SD::Expression::load (Archive &archive, const unsigned int version) |
template<class Archive > | |
void | Differentiation::SD::Expression::serialize (Archive &archive, const unsigned int version) |
Values | |
const SymEngine::Expression & | Differentiation::SD::Expression::get_expression () const |
const SymEngine::Basic & | Differentiation::SD::Expression::get_value () const |
const SymEngine::RCP< const SymEngine::Basic > & | Differentiation::SD::Expression::get_RCP () const |
Math and relational operators with (potentially) symbolic types | |
Expression & | Differentiation::SD::Expression::operator= (const Expression &rhs) |
Expression & | Differentiation::SD::Expression::operator= (Expression &&rhs) noexcept |
Expression & | Differentiation::SD::Expression::operator+= (const Expression &rhs) |
Expression & | Differentiation::SD::Expression::operator-= (const Expression &rhs) |
Expression & | Differentiation::SD::Expression::operator*= (const Expression &rhs) |
Expression & | Differentiation::SD::Expression::operator/= (const Expression &rhs) |
Math and relational operators with numeric types | |
template<typename NumberType > | |
Expression & | Differentiation::SD::Expression::operator= (const NumberType &rhs) |
Expression | Differentiation::SD::Expression::operator- () const |
template<typename NumberType > | |
Expression & | Differentiation::SD::Expression::operator+= (const NumberType &rhs) |
template<typename NumberType > | |
Expression & | Differentiation::SD::Expression::operator-= (const NumberType &rhs) |
template<typename NumberType > | |
Expression & | Differentiation::SD::Expression::operator*= (const NumberType &rhs) |
template<typename NumberType > | |
Expression & | Differentiation::SD::Expression::operator/= (const NumberType &rhs) |
Differentiation | |
Expression | Differentiation::SD::Expression::differentiate (const Expression &symbol) const |
Expression | Differentiation::SD::Expression::differentiate (const SymEngine::RCP< const SymEngine::Symbol > &symbol) const |
Expression | Differentiation::SD::Expression::differentiate (const SymEngine::RCP< const SymEngine::Basic > &symbol) const |
Dictionary-based substitution | |
Expression | Differentiation::SD::Expression::substitute (const types::substitution_map &substitution_values) const |
Expression | Differentiation::SD::Expression::substitute (const SymEngine::map_basic_basic &substitution_values) const |
Expression | Differentiation::SD::Expression::substitute (const Expression &symbol, const Expression &value) const |
template<typename NumberType > | |
Expression | Differentiation::SD::Expression::substitute (const Expression &symbol, const NumberType &value) const |
template<typename ReturnType > | |
ReturnType | Differentiation::SD::Expression::substitute_and_evaluate (const types::substitution_map &substitution_values) const |
template<typename ReturnType > | |
ReturnType | Differentiation::SD::Expression::substitute_and_evaluate (const SymEngine::map_basic_basic &substitution_values) const |
Conversion operators | |
SymEngine::Expression | Differentiation::SD::Expression::expression |
template<typename ResultType > | |
Differentiation::SD::Expression::operator ResultType () const | |
Differentiation::SD::Expression::operator const SymEngine::Expression & () const | |
Differentiation::SD::Expression::operator const SymEngine::RCP< const SymEngine::Basic > & () const | |
SymEngine::Expression & | Differentiation::SD::Expression::get_expression () |
Independent variables | |
void | Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols (const types::substitution_map &substitution_map) |
void | Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols (const SymEngine::map_basic_basic &substitution_map) |
void | Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols (const types::symbol_vector &symbols) |
void | Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols (const SymEngine::vec_basic &symbols) |
types::symbol_vector | Differentiation::SD::BatchOptimizer< ReturnType >::get_independent_symbols () const |
std::size_t | Differentiation::SD::BatchOptimizer< ReturnType >::n_independent_variables () const |
Symbol substitution | |
void | Differentiation::SD::BatchOptimizer< ReturnType >::substitute (const types::substitution_map &substitution_map) const |
void | Differentiation::SD::BatchOptimizer< ReturnType >::substitute (const SymEngine::map_basic_basic &substitution_map) const |
void | Differentiation::SD::BatchOptimizer< ReturnType >::substitute (const types::symbol_vector &symbols, const std::vector< ReturnType > &values) const |
void | Differentiation::SD::BatchOptimizer< ReturnType >::substitute (const SymEngine::vec_basic &symbols, const std::vector< ReturnType > &values) const |
bool | Differentiation::SD::BatchOptimizer< ReturnType >::values_substituted () const |
Cell iterator functions returning ranges of iterators | |
static ::ExceptionBase & | DoFHandler< dim, spacedim >::ExcGridsDoNotMatch () |
static ::ExceptionBase & | DoFHandler< dim, spacedim >::ExcInvalidBoundaryIndicator () |
static ::ExceptionBase & | DoFHandler< dim, spacedim >::ExcNewNumbersNotConsecutive (types::global_dof_index arg1) |
static ::ExceptionBase & | DoFHandler< dim, spacedim >::ExcInvalidLevel (int arg1) |
static ::ExceptionBase & | DoFHandler< dim, spacedim >::ExcFacesHaveNoLevel () |
static ::ExceptionBase & | DoFHandler< dim, spacedim >::ExcEmptyLevel (int arg1) |
Exceptions | |
static ::ExceptionBase & | DoFTools::ExcFiniteElementsDontMatch () |
static ::ExceptionBase & | DoFTools::ExcGridNotCoarser () |
static ::ExceptionBase & | DoFTools::ExcGridsDontMatch () |
static ::ExceptionBase & | DoFTools::ExcNoFESelected () |
static ::ExceptionBase & | DoFTools::ExcInvalidBoundaryIndicator () |
Support points and interpolation | |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcShapeFunctionNotPrimitive (int arg1) |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcFENotPrimitive () |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcUnitShapeValuesDoNotExist () |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcFEHasNoSupportPoints () |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcEmbeddingVoid () |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcProjectionVoid () |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcWrongInterfaceMatrixSize (int arg1, int arg2) |
static ::ExceptionBase & | FiniteElement< dim, spacedim >::ExcInterpolationNotImplemented () |
Functions to support hp | |
static ::ExceptionBase & | FE_Q_Base< PolynomialType, dim, spacedim >::ExcFEQCannotHaveDegree0 () |
Functions which should be in DoFTools | |
static ::ExceptionBase & | FETools::ExcInvalidFEName (std::string arg1) |
static ::ExceptionBase & | FETools::ExcInvalidFEDimension (char arg1, int arg2) |
static ::ExceptionBase & | FETools::ExcInvalidFE () |
static ::ExceptionBase & | FETools::ExcFENotPrimitive () |
static ::ExceptionBase & | FETools::ExcTriangulationMismatch () |
static ::ExceptionBase & | FETools::ExcHangingNodesNotAllowed () |
static ::ExceptionBase & | FETools::ExcGridNotRefinedAtLeastOnce () |
static ::ExceptionBase & | FETools::ExcMatrixDimensionMismatch (int arg1, int arg2, int arg3, int arg4) |
static ::ExceptionBase & | FETools::ExcLeastSquaresError (double arg1) |
static ::ExceptionBase & | FETools::ExcNotGreaterThan (int arg1, int arg2) |
Access to the raw data | |
static ::ExceptionBase & | FEValuesBase< dim, spacedim >::ExcAccessToUninitializedField (std::string arg1) |
static ::ExceptionBase & | FEValuesBase< dim, spacedim >::ExcFEDontMatch () |
static ::ExceptionBase & | FEValuesBase< dim, spacedim >::ExcShapeFunctionNotPrimitive (int arg1) |
static ::ExceptionBase & | FEValuesBase< dim, spacedim >::ExcFENotPrimitive () |
Exceptions | |
static ::ExceptionBase & | Mapping< dim, spacedim >::ExcInvalidData () |
static ::ExceptionBase & | Mapping< dim, spacedim >::ExcTransformationFailed () |
static ::ExceptionBase & | Mapping< dim, spacedim >::ExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3) |
Exceptions | |
static ::ExceptionBase & | GridGenerator::ExcInvalidRadii () |
static ::ExceptionBase & | GridGenerator::ExcInvalidRepetitions (int arg1) |
static ::ExceptionBase & | GridGenerator::ExcInvalidRepetitionsDimension (int arg1) |
static ::ExceptionBase & | GridGenerator::ExcInvalidInputOrientation () |
Exceptions | |
static ::ExceptionBase & | GridTools::ExcInvalidNumberOfPartitions (int arg1) |
static ::ExceptionBase & | GridTools::ExcNonExistentSubdomain (int arg1) |
static ::ExceptionBase & | GridTools::ExcTriangulationHasBeenRefined () |
static ::ExceptionBase & | GridTools::ExcScalingFactorNotPositive (double arg1) |
template<int N> | |
static ::ExceptionBase & | GridTools::ExcPointNotFoundInCoarseGrid (Point< N > arg1) |
template<int N> | |
static ::ExceptionBase & | GridTools::ExcPointNotFound (Point< N > arg1) |
static ::ExceptionBase & | GridTools::ExcVertexNotUsed (unsigned int arg1) |
Exceptions | |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcInvalidLevel (int arg1, int arg2) |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcTriangulationNotEmpty (int arg1, int arg2) |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcGridReadError () |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcFacesHaveNoLevel () |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcEmptyLevel (int arg1) |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcNonOrientableTriangulation () |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcBoundaryIdNotFound (types::boundary_id arg1) |
static ::ExceptionBase & | Triangulation< dim, spacedim >::ExcInconsistentCoarseningFlags () |
Exceptions | |
static ::ExceptionBase & | TriaRawIterator< Accessor >::ExcDereferenceInvalidCell (Accessor arg1) |
static ::ExceptionBase & | TriaRawIterator< Accessor >::ExcDereferenceInvalidObject (Accessor arg1) |
static ::ExceptionBase & | TriaRawIterator< Accessor >::ExcAdvanceInvalidObject () |
static ::ExceptionBase & | TriaRawIterator< Accessor >::ExcInvalidComparison () |
Advancement of iterators | |
static ::ExceptionBase & | TriaIterator< Accessor >::ExcAssignmentOfUnusedObject () |
Advancement of iterators | |
static ::ExceptionBase & | TriaActiveIterator< Accessor >::ExcAssignmentOfInactiveObject () |
Cell iterator functions returning ranges of iterators | |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcNoFESelected () |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcGridsDoNotMatch () |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcInvalidBoundaryIndicator () |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcMatrixHasWrongSize (int arg1) |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcFunctionNotUseful () |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcNewNumbersNotConsecutive (types::global_dof_index arg1) |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcInvalidFEIndex (int arg1, int arg2) |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcInvalidLevel (int arg1) |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcFacesHaveNoLevel () |
static ::ExceptionBase & | hp::DoFHandler< dim, spacedim >::ExcEmptyLevel (int arg1) |
Exceptions | |
static ::ExceptionBase & | hp::FECollection< dim, spacedim >::ExcNoFiniteElements () |
static ::ExceptionBase & | hp::FECollection< dim, spacedim >::ExcNoDominatedFiniteElementAmongstChildren () |
Input/Output | |
SmartPointer< const BlockSparsityPattern, BlockSparseMatrix< number > > | BlockSparseMatrix< number >::sparsity_pattern |
static ::ExceptionBase & | BlockSparseMatrix< number >::ExcBlockDimensionMismatch () |
void | BlockSparseMatrix< number >::print_formatted (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1.) const |
Multiplications | |
static ::ExceptionBase & | FullMatrix< number >::ExcEmptyMatrix () |
static ::ExceptionBase & | FullMatrix< number >::ExcNotRegular (number arg1) |
static ::ExceptionBase & | FullMatrix< number >::ExcInvalidDestination (size_type arg1, size_type arg2, size_type arg3) |
static ::ExceptionBase & | FullMatrix< number >::ExcSourceEqualsDestination () |
static ::ExceptionBase & | FullMatrix< number >::ExcMatrixNotPositiveDefinite () |
2: Implementation of VectorSpaceVector | |
static ::ExceptionBase & | LinearAlgebra::distributed::BlockVector< Number >::ExcVectorTypeNotCompatible () |
static ::ExceptionBase & | LinearAlgebra::distributed::BlockVector< Number >::ExcIteratorRangeDoesNotMatchVectorSize () |
6: Mixed stuff | |
static ::ExceptionBase & | LinearAlgebra::distributed::Vector< Number, MemorySpace >::ExcVectorTypeNotCompatible () |
static ::ExceptionBase & | LinearAlgebra::distributed::Vector< Number, MemorySpace >::ExcNotAllowedForCuda () |
static ::ExceptionBase & | LinearAlgebra::distributed::Vector< Number, MemorySpace >::ExcNonMatchingElements (Number arg1, Number arg2, unsigned int arg3) |
static ::ExceptionBase & | LinearAlgebra::distributed::Vector< Number, MemorySpace >::ExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4) |
Construction and setup | |
Constructors, destructor, functions initializing, copying and filling an object. | |
SparsityPattern::SparsityPattern () | |
SparsityPattern::SparsityPattern (const SparsityPattern &) | |
SparsityPattern::SparsityPattern (const size_type m, const size_type n, const unsigned int max_per_row) | |
SparsityPattern::SparsityPattern (const size_type m, const size_type n, const std::vector< unsigned int > &row_lengths) | |
SparsityPattern::SparsityPattern (const size_type m, const unsigned int max_per_row) | |
SparsityPattern::SparsityPattern (const size_type m, const std::vector< unsigned int > &row_lengths) | |
SparsityPattern::SparsityPattern (const SparsityPattern &original, const unsigned int max_per_row, const size_type extra_off_diagonals) | |
SparsityPattern::~SparsityPattern () override=default | |
SparsityPattern & | SparsityPattern::operator= (const SparsityPattern &) |
virtual void | SparsityPattern::reinit (const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths) override |
void | SparsityPattern::compress () |
template<typename ForwardIterator > | |
void | SparsityPattern::copy_from (const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end) |
void | SparsityPattern::copy_from (const DynamicSparsityPattern &dsp) |
void | SparsityPattern::copy_from (const SparsityPattern &sp) |
template<typename number > | |
void | SparsityPattern::copy_from (const FullMatrix< number > &matrix) |
template<typename ForwardIterator > | |
void | SparsityPattern::add_entries (const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false) |
Querying information | |
bool | SparsityPattern::operator== (const SparsityPattern &) const |
bool | SparsityPattern::stores_only_added_elements () const |
std::size_t | SparsityPattern::memory_consumption () const |
Accessing entries | |
size_type | SparsityPattern::operator() (const size_type i, const size_type j) const |
Partitioners | |
static ::ExceptionBase & | TrilinosWrappers::PreconditionBase::ExcNonMatchingMaps (std::string arg1) |
Input/Output | |
static ::ExceptionBase & | TrilinosWrappers::SparseMatrix::ExcTrilinosError (int arg1) |
static ::ExceptionBase & | TrilinosWrappers::SparseMatrix::ExcInvalidIndex (size_type arg1, size_type arg2) |
static ::ExceptionBase & | TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination () |
static ::ExceptionBase & | TrilinosWrappers::SparseMatrix::ExcMatrixNotCompressed () |
static ::ExceptionBase & | TrilinosWrappers::SparseMatrix::ExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4) |
static ::ExceptionBase & | TrilinosWrappers::SparseMatrix::ExcAccessToNonPresentElement (size_type arg1, size_type arg2) |
4: Mixed stuff | |
static ::ExceptionBase & | TrilinosWrappers::MPI::Vector::ExcDifferentParallelPartitioning () |
static ::ExceptionBase & | TrilinosWrappers::MPI::Vector::ExcTrilinosError (int arg1) |
static ::ExceptionBase & | TrilinosWrappers::MPI::Vector::ExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4) |
This module contains classes that are used in the exception mechanism of deal.II.
Exceptions are used in two different ways:
Static assertions: These are checks that are only enabled in debug mode, not in release (or optimized, production) mode. In deal.II, static assertions are typically used to check that parameters to functions satisfy certain properties, that internal data structures are consistent, and similar assertions. For example, static assertions are used to make sure that two vectors that are added together have the same number of components – everything else would not make any sense anyway.
Such checks are performed by the Assert
macro in several thousand places within the library. Also, several tutorial programs starting with step-5 show how to do this.
If a static assertion is violated, the exception mechanism generates an exception of a type that indicates what exactly goes wrong, displays appropriate information including the exact location where the problem was detected, and then aborts the program – if you try to add two vectors of different length, there is nothing that can be done within the program to cope with the situation, you have to go fix the program code instead. There is generally not even a reason to throw
an exception object using the usual C++ exception mechanism because there is nothing a function higher up could do in such cases to rectify the situation and deal with it in a useful way – it's not that the program received bad data; the program is just buggy, and one can not intelligently work around that.
(It is sometimes useful to change the behavior of the Assert
macro from aborting a program to throwing exceptions. On the other hand, exceptions are not allowed to propagate out of destructors of classes. For this purpose, there is a variant of the macro, called AssertNothrow
, that can be used in destructors. These use cases are discussed further down below on this page.)
Dynamic assertions: These are used to check conditions that depend on external things that may be different from one program run to the next, such as whether an output file can be written to.
These are things that shouldn't be checked statically, because it is not guaranteed that a program for which the condition is satisfied in a debug mode run, will also have the condition satisfied in a subsequent release mode run – in other words, it is not sufficient to only check these situations in debug mode.
Rather, one has to check them every time during execution of a program. Within deal.II, this is done using the AssertThrow
macro introduced in step-9, step-13, and following tutorial programs. The macro checks a condition, and if violated throws an exception of one of the types declared in this module, using the C++ throw
mechanism. Since these are run-time exceptions, this gives the program the chance to catch the exception and, for example, write the output to a writable file instead.
The error handling mechanism in deal.II
is generally used in two ways. The first uses error checking in debug mode only and is useful for programs which are not fully tested. When the program shows no errors anymore, one may switch off error handling and get better performance by this, since checks for errors are done quite frequently in the library (a typical speed up is a factor of four!). This mode of exception generation is most useful for internal consistency checks such as range checking or checking of the validity of function arguments. Errors of this kind usually are programming errors and the program should abort with as detailed a message as possible, including location and reason for the generation of the exception.
The second mode is for error checks which should always be on, such as for I/O errors, failing memory requests and the like. It does not make much sense to turn this mode off, since this kind of errors may happen in tested and untested programs likewise. Exceptions of this kind do not terminate the program, rather they throw exceptions in the C++
manner, allowing the program to catch them and eventually do something about it. As it may be useful to have some information printed out if an exception could not be handled properly, additional information is passed along as for the first mode. The latter makes it necessary to provide a family of macros which enter this additional information into the exception class; this could in principle be done by the programmer himself each time by hand, but since the information can be obtained automatically, a macro is provided for this.
Both modes use exception classes, which need to have special features in additional to the C++
standard's std::exception
class. Such a class is declared by the following lines of code:
This declares an exception class named ExcDomain
, which has two variables as additional information (named arg1
and arg2
by default) and which outputs the given sequence (which is appended to an std::ostream
variable's name, thus the weird syntax). There are other DeclExceptionN
macros for exception classes with more or no parameters. By convention, the name of all exception classes starts with Exc...
and most of them are declared locally to the class it is to be used in (a few very frequently found ones are also declared in the StandardExceptions namespace and are available everywhere). Declaring exceptions globally is possible but pollutes the global namespace, is less readable and most of the time unnecessary.
Since exception classes are declared the same way for both modes of error checking, it is possible to use an exception declared through the DeclExceptionN(...)
macro family for both static as well as dynamic checks.
To use the exception mechanism for debug mode error checking, write lines like the following in your source code:
which by macro expansion does essentially the following (though the actual code is slightly more complicated):
That is, it issues an error only if the preprocessor variable DEBUG
is set and if the given condition (in this case n < dim
is violated).
If the exception was declared using the DeclException0 (...)
macro, i.e., without any additional parameters, its name has nonetheless to be given with parentheses: Assert (i>m, ExcSomewhat());
If the DEBUG
preprocessor directive is set, the call Assert (cond, exc);
is basically converted by the preprocessor into the following sequence:
(Note that function names and exact calling sequences may change over time, but the general principle remains the same.) I.e., if the given condition is violated, then the file and line in which the exception occurred as well as the condition itself and the call sequence of the exception object is passed to the deal_II_exceptions::internals::issue_error_noreturn() function. Additionally an object of the form given by exc
is created (this is normally an unnamed object like in ExcDomain (n, dim)
of class ExcDomain
) and transferred to this function.
__PRETTY_FUNCTION__
is a macro defined by some compilers and gives the name of the function. If another compiler is used, we try to set this function to something reasonable, if the compiler provides us with that, and "(not available)"
otherwise.
In issue_error_noreturn
, the given data is transferred into the exc
object by calling the set_fields() function; Afterwards the program is either aborted (and information about the exception is printed to deallog) or the exception is thrown. The Assert
macro does the first path (print and abort); AssertThrow
does the second (throw). This behavior is consistent with the descriptions of static and dynamic assertions earlier in this document. If it can be obtained from the operating system, the output may also contain a stacktrace to show where the error happened. Several of the Tutorial programs programs show a typical output.
If the preprocessor variable DEBUG
is not set then the Assert
macro is expanded to {}
.
Sometimes, there is no useful condition for an exception other than that the program flow should not have reached a certain point, e.g. a default
section of a switch
statement. In this case, raise the exception by the following construct:
See the step-7 and several other of the tutorial programs for a use of this construct.
As mentioned above, the program is terminated once a call to Assert
fails. However, there is one case where we do not want to do this, namely when a C++ exception is active. The usual case where this happens is that someone throws an exception through the AssertThrow
mechanism (see below) which, while the stack is unwound, leads to the destruction of other objects in stack frames above. If other objects refer to the objects being thus destroyed, some destructors raise an exception through Assert
. If we would abort the program then, we would only ever see the message that an object is being destroyed which is still referenced from somewhere, but we would never see the original exception that triggered this. (You can see it in the debugger by putting a break point on the function __throw
, but you cannot see it from the program itself.) In that case, we use a C++ standard library function to detect the presence of another active exception and do not terminate the program to allow that the thrown exception propagates to some place where its message can be displayed.
Since it is common that one failed assertion leads to a whole chain of others, we only ever print the very first message. If the program is then aborted, that is no problem. If it is not (since a C++ exception is active), only the first is displayed and a message about suppressed follow-up messages is shown.
C++ has a mechanism to indicate that something exceptional has happened: exceptions that can be triggered by throw
statements and captured by catch
clauses, see for example https://en.wikipedia.org/wiki/C%2B%2B#Exception_handling and http://www.cplusplus.com/doc/tutorial/exceptions/ .
At some fundamental level, a typical C++ exception is an object that is placed in some special place, and then the function exits the current scope (e.g., the current function) through an exceptional return path. This is often enough to tell what problem triggered the exception, but more frequently it would be nice if one had more information: for example, in which line of the code the problem happened, or what non-existent entry of a sparse matrix the code wanted to write into.
Dynamic assertions in deal.II therefore extend this mechanism a bit. Typically, one would raise an exception by code such as
and catch it using the statement
std::exception
is a standard C++
class providing basic functionality for exceptions, such as the virtual function what()
that returns some information on the exception itself. This information is useful if an exception can't be handled properly, in which case as precise a description as possible should be printed.
The problem here is that to get significant and useful information out of what()
, it is necessary to overload this function in our exception class and call the throw
operator with additional arguments to the exception class. The first thing, overloading the what
function is done using the DeclExceptionN
macros, but putting the right information, which is the same as explained above for the Assert
expansion, requires some work if one would want to write it down each time:
For this purpose, the macro AssertThrow
was invented. It does mainly the same job as does the Assert
macro, but it does not abort the program; rather, it throws an exception as shown above. The mode of usage is
The condition to be checked is incorporated into the macro in order to allow passing the violated condition as a string. The expansion of the AssertThrow
macro is not affected by the DEBUG
preprocessor variable.
There is a whole family of DeclExceptionX
macros where X
is to be replaced by the number of additional parameters (0 to 5 presently). These macros are used to declare exception classes in the following way:
The first argument denotes the name of the exception class to be created. The next arguments are the types of the parameters (in this case there two types, corresponding to the X
in DeclExceptionX
) and finally the output sequence with which you can print additional information.
The syntax of the output sequence is a bit weird but gets clearer once you see how this macro is defined (again schematically, actual function names and definitions may change over time and be different):
If declared as specified, you can later use this exception class in the following manner:
and the output if the condition fails will be
Obviously for the DeclException0(name)
macro, no types and also no output sequence is allowed.
Assert:
The AssertNothrow
macroThe default implementation of the Assert
macro, as discussed above, prints detailed information about what exactly went wrong to the screen and then aborts the program. Aborting the program is useful because it allows easily finding the place where something went wrong – including all of the information how we got to that place – by running the program in a debugger.
On the other hand, there are cases where aborting a program may be undesirable and one needs to exit in a somewhat more graceful way – even if there is really not very much one can do in these cases to still produce a meaningful result. An example is if a deal.II program is run a one module in a bigger framework of software. Think, for example, of a case where a deal.II program computed the flow field that corresponds to a set of input variables provided by some optimization routine: if the optimizer on the outside provided a negative density as input (a condition one might want to check via Assert
), then this clearly makes no sense, and the flow solver cannot produce a meaningful answer; but it should tell that to the optimizer nicely, rather than just aborting the entire process (optimizer and flow solver).
For this purpose, one can call deal_II_exceptions::disable_abort_on_exception() that switches what Assert
does from aborting the program to essentially the same as AssertThrow
does, namely using the C++ throw
mechanism to raise an exception. This exception can then be caught at a higher level – e.g., in the optimization routine that sits atop the flow solver, and that can then decide what it wants to do with the situation.
This is all nice and good, but C++ does not allow throwing exceptions inside the destructors of classes, or in a function that is currently being called from a destructor higher up in the call stack. To this end, there is a separate macro, AssertNothrow
, that can be used in destructors: It acts just like Assert
usually does – in particular, it only checks the condition in debug mode – but it is immune to the effect of deal_II_exceptions::disable_abort_on_exception(): It will only ever abort the program, and never throw an exception.
#define DeclException0 | ( | Exception0 | ) |
Declare an exception class derived from ExceptionBase without parameters.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 473 of file exceptions.h.
#define DeclExceptionMsg | ( | Exception, | |
defaulttext | |||
) |
Declare an exception class derived from ExceptionBase that can take one runtime argument, but if none is given in the place where you want to throw the exception, it simply reverts to the default text provided when declaring the exception class through this macro.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included. Definition at line 496 of file exceptions.h.
#define DeclException1 | ( | Exception1, | |
type1, | |||
outsequence | |||
) |
Declare an exception class derived from ExceptionBase with one additional parameter.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included. Definition at line 518 of file exceptions.h.
#define DeclException2 | ( | Exception2, | |
type1, | |||
type2, | |||
outsequence | |||
) |
Declare an exception class derived from ExceptionBase with two additional parameters.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included. Definition at line 541 of file exceptions.h.
#define DeclException3 | ( | Exception3, | |
type1, | |||
type2, | |||
type3, | |||
outsequence | |||
) |
Declare an exception class derived from ExceptionBase with three additional parameters.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included. Definition at line 564 of file exceptions.h.
#define DeclException4 | ( | Exception4, | |
type1, | |||
type2, | |||
type3, | |||
type4, | |||
outsequence | |||
) |
Declare an exception class derived from ExceptionBase with four additional parameters.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included. Definition at line 587 of file exceptions.h.
#define DeclException5 | ( | Exception5, | |
type1, | |||
type2, | |||
type3, | |||
type4, | |||
type5, | |||
outsequence | |||
) |
Declare an exception class derived from ExceptionBase with five additional parameters.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included. Definition at line 613 of file exceptions.h.
#define Assert | ( | cond, | |
exc | |||
) |
A macro that serves as the main routine in the exception mechanism for debug mode error checking. It asserts that a certain condition is fulfilled, otherwise issues an error and aborts the program.
A more detailed description can be found in the Exceptions and assertions module. It is first used in step-5 and step-6. See also the ExceptionBase
class for more information.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1419 of file exceptions.h.
#define AssertNothrow | ( | cond, | |
exc | |||
) |
A variant of the Assert
macro above that exhibits the same runtime behavior as long as disable_abort_on_exception was not called.
However, if disable_abort_on_exception was called, this macro merely prints the exception that would be thrown to deallog and continues normally without throwing an exception.
A more detailed description can be found in the Exceptions and assertions module, in the discussion about the corner case at the bottom of the page.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1483 of file exceptions.h.
#define AssertThrow | ( | cond, | |
exc | |||
) |
A macro that serves as the main routine in the exception mechanism for dynamic error checking. It asserts that a certain condition is fulfilled, otherwise throws an exception via the C++ throw
mechanism. This exception can be caught via a catch
clause, as is shown in step-6 and all following tutorial programs.
A more detailed description can be found in the Exceptions and assertions module. It is first used in step-9 and step-13. See also the ExceptionBase
class for more information.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1531 of file exceptions.h.
#define AssertDimension | ( | dim1, | |
dim2 | |||
) |
Special assertion for dimension mismatch.
Since this is used very often and always repeats the arguments, we introduce this special assertion for ExcDimensionMismatch in order to keep the user codes shorter.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1579 of file exceptions.h.
#define AssertVectorVectorDimension | ( | VEC, | |
DIM1, | |||
DIM2 | |||
) |
An assertion that tests whether vec
has size dim1
, and each entry of the vector is itself an array that has the size dim2
.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1607 of file exceptions.h.
#define AssertIndexRange | ( | index, | |
range | |||
) |
An assertion that tests that a given index is within the half-open range [0,range)
. It throws an exception object ExcIndexRange(index,0,range)
if the assertion fails.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1649 of file exceptions.h.
#define AssertIsFinite | ( | number | ) |
An assertion that checks whether a number is finite or not. We explicitly cast the number to std::complex to match the signature of the exception (see there for an explanation of why we use std::complex at all) and to satisfy the fact that std::complex has no implicit conversions.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1681 of file exceptions.h.
#define AssertThrowMPI | ( | error_code | ) | AssertThrow(error_code == MPI_SUCCESS, ::ExcMPI(error_code)) |
An assertion that checks whether or not an error code returned by an MPI function is equal to MPI_SUCCESS
. If the check fails then an exception of type ExcMPI is thrown with the given error code as an argument.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1707 of file exceptions.h.
#define AssertCuda | ( | error_code | ) |
An assertion that checks that the error code produced by calling a CUDA routine is equal to cudaSuccess.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1734 of file exceptions.h.
#define AssertNothrowCuda | ( | error_code | ) |
The non-throwing equivalent of AssertCuda.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1762 of file exceptions.h.
#define AssertCudaKernel | ( | ) |
An assertion that checks that the kernel was launched and executed successfully.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1791 of file exceptions.h.
#define AssertCusparse | ( | error_code | ) |
An assertion that checks that the error code produced by calling a cuSPARSE routine is equal to CUSPARSE_STATUS_SUCCESS.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1822 of file exceptions.h.
#define AssertNothrowCusparse | ( | error_code | ) |
The non-throwing equivalent of AssertCusparse.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1853 of file exceptions.h.
#define AssertCusolver | ( | error_code | ) |
An assertion that checks that the error code produced by calling a cuSOLVER routine is equal to CUSOLVER_STATUS_SUCCESS.
#undef
this macro, as well as all other macros defined by deal.II that are not prefixed with either DEAL
or deal
, by including the header deal.II/base/undefine_macros.h
after all other deal.II headers have been included.Definition at line 1885 of file exceptions.h.
using DataOutBase::EpsFlags::ColorFunction = RgbValues (*)(const double value, const double min_value, const double max_value) |
Definition of a function pointer type taking a value and returning a triple of color values in RGB values.
Besides the actual value by which the color is to be computed, min and max values of the data to be colorized are given as well.
Definition at line 953 of file data_out_base.h.
|
private |
Definition at line 1482 of file data_out_base.h.
Exception for when an entry in an XML parameter file does not match the provided pattern. The arguments are, in order, the entry value, entry name, and a description of the pattern.
Definition at line 1712 of file parameter_handler.h.
using Patterns::Tools::Convert< Tensor< rank, dim, Number > >::T = Tensor<rank, dim, Number> |
Definition at line 1935 of file patterns.h.
using Patterns::Tools::Convert< Point< dim, Number > >::T = Point<dim, Number> |
Definition at line 2002 of file patterns.h.
using Patterns::Tools::Convert< std::unique_ptr< FunctionParser< dim > > >::T = std::unique_ptr<FunctionParser<dim> > |
Definition at line 2032 of file patterns.h.
using Patterns::Tools::Convert< ComponentMask >::T = ComponentMask |
Definition at line 2098 of file patterns.h.
using Patterns::Tools::Convert< std::complex< Number > >::T = std::complex<Number> |
Definition at line 2132 of file patterns.h.
using Patterns::Tools::Convert< std::string >::T = std::string |
Definition at line 2196 of file patterns.h.
using Patterns::Tools::Convert< std::pair< Key, Value > >::T = std::pair<Key, Value> |
Definition at line 2227 of file patterns.h.
using Patterns::Tools::Convert< std::tuple< Args... > >::T = std::tuple<Args...> |
Definition at line 2266 of file patterns.h.
|
private |
The data type used in counter_map.
Definition at line 222 of file subscriptor.h.
|
private |
The iterator type used in counter_map.
Definition at line 227 of file subscriptor.h.
|
private |
A map type used to indicate which dependent variable is associated with which entry in the output vector.
Definition at line 1998 of file symengine_optimizer.h.
using BlockVectorBase< VectorType >::BlockType = VectorType |
Typedef the type of the underlying vector.
Definition at line 456 of file block_vector_base.h.
using BlockVectorBase< VectorType >::value_type = typename BlockType::value_type |
Definition at line 468 of file block_vector_base.h.
using BlockVectorBase< VectorType >::pointer = value_type * |
Definition at line 469 of file block_vector_base.h.
using BlockVectorBase< VectorType >::const_pointer = const value_type * |
Definition at line 470 of file block_vector_base.h.
using BlockVectorBase< VectorType >::iterator = ::internal::BlockVectorIterators::Iterator<BlockVectorBase, false> |
Definition at line 472 of file block_vector_base.h.
using BlockVectorBase< VectorType >::const_iterator = ::internal::BlockVectorIterators::Iterator<BlockVectorBase, true> |
Definition at line 474 of file block_vector_base.h.
using BlockVectorBase< VectorType >::reference = typename BlockType::reference |
Definition at line 475 of file block_vector_base.h.
using BlockVectorBase< VectorType >::const_reference = typename BlockType::const_reference |
Definition at line 476 of file block_vector_base.h.
using BlockVectorBase< VectorType >::size_type = types::global_dof_index |
Definition at line 477 of file block_vector_base.h.
using BlockVectorBase< VectorType >::real_type = typename BlockType::real_type |
Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. If the template argument of this class is a real data type, then real_type equals the template argument. If the template argument is a std::complex type then real_type equals the type underlying the complex numbers.
This alias is used to represent the return type of norms.
Definition at line 488 of file block_vector_base.h.
|
private |
Define number type of matrix.
Definition at line 390 of file precondition_block.h.
using PreconditionBlockJacobi< MatrixType, inverse_type >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 396 of file precondition_block.h.
using PreconditionBlockSOR< MatrixType, inverse_type >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 667 of file precondition_block.h.
using PreconditionBlockSOR< MatrixType, inverse_type >::number = typename MatrixType::value_type |
Define number type of matrix.
Definition at line 677 of file precondition_block.h.
using PreconditionBlockSSOR< MatrixType, inverse_type >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 832 of file precondition_block.h.
using PreconditionBlockSSOR< MatrixType, inverse_type >::number = typename MatrixType::value_type |
Define number type of matrix.
Definition at line 837 of file precondition_block.h.
Declare type for container size.
Definition at line 871 of file sparsity_pattern.h.
Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Definition at line 877 of file sparsity_pattern.h.
Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Since the iterator does not allow to modify the sparsity pattern, this type is the same as that for const_iterator
.
Definition at line 886 of file sparsity_pattern.h.
using GrowingVectorMemory< VectorType >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 326 of file vector_memory.h.
|
private |
A type that describes this entries of an array that represents the vectors stored by this object. The first component of the pair is be a flag telling whether the vector is used, the second a pointer to the vector itself.
Definition at line 403 of file vector_memory.h.
Enum denoting the possibilities whether the scaling should be done such that the given size
equals the width or the height of the resulting picture.
Enumerator | |
---|---|
width | Scale to given width. |
height | Scale to given height. |
Definition at line 820 of file data_out_base.h.
A data type providing the different possible zlib compression levels. These map directly to constants defined by zlib.
Definition at line 1133 of file data_out_base.h.
Declare what a multiple entry is: a variant entry (in curly braces {
, }
) or an array (in double curly braces {{
, }}
).
Enumerator | |
---|---|
variant | A variant entry. |
array | An array entry. |
Definition at line 2185 of file parameter_handler.h.
Files can be used for input or output. This can be specified in the constructor by choosing the flag type
.
Enumerator | |
---|---|
input | Open for input. |
output | Open for output. |
Definition at line 1093 of file patterns.h.
|
strong |
An enumeration to distinguish between different optimization methods that can be used by SymEngine to more rapidly evaluate complex symbolic expressions.
Definition at line 99 of file symengine_optimizer.h.
|
strong |
An enumeration to specify which special techniques, over and above those used with the chosen OptimizerType, to be applied to the set of expressions that are to be optimized.
Definition at line 145 of file symengine_optimizer.h.
|
static |
An entry with this name does not exist in the AnyData object.
|
static |
The requested type and the stored type are different.
|
static |
Exception indicating that a function expected a vector to have a certain name, but we store a different name in that position.
|
static |
An entry with this name does not exist in the internal boost::any map.
|
static |
An entry with this name does not exist in the internal boost::any map.
|
static |
The requested type and the stored type are different.
|
static |
Exception
DataOutBase::Patch< dim, spacedim >::Patch |
Default constructor. Sets n_subdivisions to one, points_are_available to false, and patch_index to no_neighbor.
Definition at line 1611 of file data_out_base.cc.
bool DataOutBase::Patch< dim, spacedim >::operator== | ( | const Patch< dim, spacedim > & | patch | ) | const |
Compare the present patch for equality with another one. This is used in a few of the automated tests in our testsuite.
Definition at line 1629 of file data_out_base.cc.
std::size_t DataOutBase::Patch< dim, spacedim >::memory_consumption |
Return an estimate for the memory consumption, in bytes, of this object. This is not exact (but will usually be close) because calculating the memory usage of trees (e.g., std::map
) is difficult.
Definition at line 1668 of file data_out_base.cc.
void DataOutBase::Patch< dim, spacedim >::swap | ( | Patch< dim, spacedim > & | other_patch | ) |
Swap the current object's contents with those of the given argument.
Definition at line 1684 of file data_out_base.cc.
|
static |
Exception
DataOutBase::Patch< 0, spacedim >::Patch |
Default constructor. Sets points_are_available to false, and patch_index to no_neighbor.
Definition at line 1712 of file data_out_base.cc.
bool DataOutBase::Patch< 0, spacedim >::operator== | ( | const Patch< 0, spacedim > & | patch | ) | const |
Compare the present patch for equality with another one. This is used in a few of the automated tests in our testsuite.
Definition at line 1723 of file data_out_base.cc.
std::size_t DataOutBase::Patch< 0, spacedim >::memory_consumption |
Return an estimate for the memory consumption, in bytes, of this object. This is not exact (but will usually be close) because calculating the memory usage of trees (e.g., std::map
) is difficult.
Definition at line 1757 of file data_out_base.cc.
void DataOutBase::Patch< 0, spacedim >::swap | ( | Patch< 0, spacedim > & | other_patch | ) |
Swap the current object's contents with those of the given argument.
Definition at line 1768 of file data_out_base.cc.
|
static |
Exception
|
static |
Declare all flags with name and type as offered by this class, for use in input files.
This method does nothing, but child classes may override this method to add fields to prm
.
Definition at line 571 of file data_out_base.h.
void DataOutBase::OutputFlagsBase< FlagsType >::parse_parameters | ( | const ParameterHandler & | prm | ) |
Read the parameters declared in declare_parameters() and set the flags for this output format accordingly.
This method does nothing, but child classes may override this method to add fields to prm
.
Definition at line 577 of file data_out_base.h.
std::size_t DataOutBase::OutputFlagsBase< FlagsType >::memory_consumption |
Return an estimate for the memory consumption, in bytes, of this object. This is not exact (but will usually be close) because calculating the memory usage of trees (e.g., std::map
) is difficult.
Definition at line 583 of file data_out_base.h.
DataOutBase::DXFlags::DXFlags | ( | const bool | write_neighbors = false , |
const bool | int_binary = false , |
||
const bool | coordinates_binary = false , |
||
const bool | data_binary = false |
||
) |
Constructor.
Definition at line 1870 of file data_out_base.cc.
|
static |
Declare all flags with name and type as offered by this class, for use in input files.
Definition at line 1883 of file data_out_base.cc.
void DataOutBase::DXFlags::parse_parameters | ( | const ParameterHandler & | prm | ) |
Read the parameters declared in declare_parameters() and set the flags for this output format accordingly.
The flags thus obtained overwrite all previous contents of this object.
Definition at line 1914 of file data_out_base.cc.
DataOutBase::UcdFlags::UcdFlags | ( | const bool | write_preamble = false | ) |
Constructor.
Definition at line 1778 of file data_out_base.cc.
|
static |
Declare all flags with name and type as offered by this class, for use in input files.
Definition at line 1923 of file data_out_base.cc.
void DataOutBase::UcdFlags::parse_parameters | ( | const ParameterHandler & | prm | ) |
Read the parameters declared in declare_parameters() and set the flags for this output format accordingly.
The flags thus obtained overwrite all previous contents of this object.
Definition at line 1937 of file data_out_base.cc.