class SymmetricTensor;
template <int dim, typename Number>
-SymmetricTensor<2, dim, Number>
-unit_symmetric_tensor();
+DEAL_II_CONSTEXPR SymmetricTensor<2, dim, Number>
+ unit_symmetric_tensor();
template <int dim, typename Number>
-SymmetricTensor<4, dim, Number>
-deviator_tensor();
+DEAL_II_CONSTEXPR SymmetricTensor<4, dim, Number>
+ deviator_tensor();
template <int dim, typename Number>
-SymmetricTensor<4, dim, Number>
-identity_tensor();
+DEAL_II_CONSTEXPR SymmetricTensor<4, dim, Number>
+ identity_tensor();
template <int dim, typename Number>
-SymmetricTensor<2, dim, Number>
+constexpr SymmetricTensor<2, dim, Number>
invert(const SymmetricTensor<2, dim, Number> &);
template <int dim, typename Number>
-SymmetricTensor<4, dim, Number>
+constexpr SymmetricTensor<4, dim, Number>
invert(const SymmetricTensor<4, dim, Number> &);
template <int dim2, typename Number>
-Number
-trace(const SymmetricTensor<2, dim2, Number> &);
+DEAL_II_CONSTEXPR Number
+ trace(const SymmetricTensor<2, dim2, Number> &);
template <int dim, typename Number>
-SymmetricTensor<2, dim, Number>
-deviator(const SymmetricTensor<2, dim, Number> &);
+DEAL_II_CONSTEXPR SymmetricTensor<2, dim, Number>
+ deviator(const SymmetricTensor<2, dim, Number> &);
template <int dim, typename Number>
-Number
-determinant(const SymmetricTensor<2, dim, Number> &);
+DEAL_II_CONSTEXPR Number
+ determinant(const SymmetricTensor<2, dim, Number> &);
* put at position <tt>position</tt>. The remaining indices remain in
* invalid state.
*/
- inline TableIndices<2>
+ DEAL_II_CONSTEXPR inline TableIndices<2>
merge(const TableIndices<2> &previous_indices,
const unsigned int new_index,
const unsigned int position)
* put at position <tt>position</tt>. The remaining indices remain in
* invalid state.
*/
- inline TableIndices<4>
+ DEAL_II_CONSTEXPR inline TableIndices<4>
merge(const TableIndices<4> &previous_indices,
const unsigned int new_index,
const unsigned int position)
previous_indices[1],
previous_indices[2],
new_index};
+ default:
+ Assert(false, ExcInternalError());
+ return {};
}
- Assert(false, ExcInternalError());
- return {};
}
* @internal
*
* Class that acts as accessor to elements of type SymmetricTensor. The
- * template parameter <tt>C</tt> may be either true or false, and
+ * template parameter <tt>constness</tt> may be either true or false, and
* indicates whether the objects worked on are constant or not (i.e. write
* access is only allowed if the value is false).
*
* Since with <tt>N</tt> indices, the effect of applying
- * <tt>operator[]</tt> is getting access to something we <tt>N-1</tt>
+ * <tt>operator[]</tt> is getting access to something with <tt>N-1</tt>
* indices, we have to implement these accessor classes recursively, with
* stopping when we have only one index left. For the latter case, a
* specialization of this class is declared below, where calling
* This guarantees that the accessor objects go out of scope earlier
* than the mother object, avoid problems with data consistency.
*/
- Accessor(tensor_type &tensor, const TableIndices<rank> &previous_indices);
+ constexpr Accessor(tensor_type & tensor,
+ const TableIndices<rank> &previous_indices);
/**
* Copy constructor.
*/
- Accessor(const Accessor &) = default;
+ constexpr Accessor(const Accessor &) = default;
public:
/**
* Index operator.
*/
- Accessor<rank, dim, constness, P - 1, Number>
- operator[](const unsigned int i);
+ DEAL_II_CONSTEXPR Accessor<rank, dim, constness, P - 1, Number>
+ operator[](const unsigned int i);
/**
* Index operator.
*/
- Accessor<rank, dim, constness, P - 1, Number>
+ constexpr Accessor<rank, dim, constness, P - 1, Number>
operator[](const unsigned int i) const;
private:
* This guarantees that the accessor objects go out of scope earlier
* than the mother object, avoid problems with data consistency.
*/
- Accessor(tensor_type &tensor, const TableIndices<rank> &previous_indices);
+ constexpr Accessor(tensor_type & tensor,
+ const TableIndices<rank> &previous_indices);
/**
* Copy constructor.
*/
- Accessor(const Accessor &) = default;
+ constexpr Accessor(const Accessor &) = default;
public:
/**
* Index operator.
*/
- reference operator[](const unsigned int);
+ DEAL_II_CONSTEXPR reference operator[](const unsigned int);
/**
* Index operator.
*/
- reference operator[](const unsigned int) const;
+ constexpr reference operator[](const unsigned int) const;
private:
/**
/**
* Default constructor. Creates a tensor with all entries equal to zero.
*/
- SymmetricTensor();
+ constexpr SymmetricTensor() = default;
/**
* Constructor. Generate a symmetric tensor from a general one. Assumes that
* practice to check before calling <tt>symmetrize</tt>.
*/
template <typename OtherNumber>
- explicit SymmetricTensor(const Tensor<2, dim, OtherNumber> &t);
+ DEAL_II_CONSTEXPR explicit SymmetricTensor(
+ const Tensor<2, dim, OtherNumber> &t);
/**
* A constructor that creates a symmetric tensor from an array holding its
* the object from the internal namespace is to work around bugs in some
* older compilers.
*/
+ DEAL_II_CONSTEXPR
SymmetricTensor(const Number (&array)[n_independent_components]);
/**
* Number.
*/
template <typename OtherNumber>
- explicit SymmetricTensor(
+ constexpr explicit SymmetricTensor(
const SymmetricTensor<rank_, dim, OtherNumber> &initializer);
/**
* Return a pointer to the first element of the underlying storage.
*/
- Number *
- begin_raw();
+ DEAL_II_CONSTEXPR Number *
+ begin_raw();
/**
* Return a const pointer to the first element of the underlying storage.
*/
- const Number *
+ constexpr const Number *
begin_raw() const;
/**
* Return a pointer to the element past the end of the underlying storage.
*/
- Number *
- end_raw();
+ DEAL_II_CONSTEXPR Number *
+ end_raw();
/**
* Return a const pointer to the element past the end of the underlying
* storage.
*/
- const Number *
+ constexpr const Number *
end_raw() const;
/**
* @p Number.
*/
template <typename OtherNumber>
- SymmetricTensor &
- operator=(const SymmetricTensor<rank_, dim, OtherNumber> &rhs);
+ DEAL_II_CONSTEXPR SymmetricTensor &
+ operator=(const SymmetricTensor<rank_, dim, OtherNumber> &rhs);
/**
* This operator assigns a scalar to a tensor. To avoid confusion with what
* value allowed for <tt>d</tt>, allowing the intuitive notation
* <tt>t=0</tt> to reset all elements of the tensor to zero.
*/
- SymmetricTensor &
- operator=(const Number &d);
+ DEAL_II_CONSTEXPR SymmetricTensor &
+ operator=(const Number &d);
/**
* Convert the present symmetric tensor into a full tensor with the same
* elements, but using the different storage scheme of full tensors.
*/
- operator Tensor<rank_, dim, Number>() const;
+ constexpr operator Tensor<rank_, dim, Number>() const;
/**
* Test for equality of two tensors.
*/
- bool
+ constexpr bool
operator==(const SymmetricTensor &) const;
/**
* Test for inequality of two tensors.
*/
- bool
+ constexpr bool
operator!=(const SymmetricTensor &) const;
/**
* Add another tensor.
*/
template <typename OtherNumber>
- SymmetricTensor &
- operator+=(const SymmetricTensor<rank_, dim, OtherNumber> &);
+ DEAL_II_CONSTEXPR SymmetricTensor &
+ operator+=(const SymmetricTensor<rank_, dim, OtherNumber> &);
/**
* Subtract another tensor.
*/
template <typename OtherNumber>
- SymmetricTensor &
- operator-=(const SymmetricTensor<rank_, dim, OtherNumber> &);
+ DEAL_II_CONSTEXPR SymmetricTensor &
+ operator-=(const SymmetricTensor<rank_, dim, OtherNumber> &);
/**
* Scale the tensor by <tt>factor</tt>, i.e. multiply all components by
* <tt>factor</tt>.
*/
template <typename OtherNumber>
- SymmetricTensor &
- operator*=(const OtherNumber &factor);
+ DEAL_II_CONSTEXPR SymmetricTensor &
+ operator*=(const OtherNumber &factor);
/**
* Scale the tensor by <tt>1/factor</tt>.
*/
template <typename OtherNumber>
- SymmetricTensor &
- operator/=(const OtherNumber &factor);
+ DEAL_II_CONSTEXPR SymmetricTensor &
+ operator/=(const OtherNumber &factor);
/**
* Unary minus operator. Negate all entries of a tensor.
*/
- SymmetricTensor
- operator-() const;
+ DEAL_II_CONSTEXPR SymmetricTensor
+ operator-() const;
/**
* Product between the present symmetric tensor and a tensor of rank 2. For
* they write it into the first argument to the function.
*/
template <typename OtherNumber>
- typename internal::SymmetricTensorAccessors::
+ DEAL_II_CONSTEXPR typename internal::SymmetricTensorAccessors::
double_contraction_result<rank_, 2, dim, Number, OtherNumber>::type
operator*(const SymmetricTensor<2, dim, OtherNumber> &s) const;
* symmetric tensor given as argument.
*/
template <typename OtherNumber>
- typename internal::SymmetricTensorAccessors::
+ DEAL_II_CONSTEXPR typename internal::SymmetricTensorAccessors::
double_contraction_result<rank_, 4, dim, Number, OtherNumber>::type
operator*(const SymmetricTensor<4, dim, OtherNumber> &s) const;
/**
* Return a read-write reference to the indicated element.
*/
- Number &
- operator()(const TableIndices<rank_> &indices);
+ DEAL_II_CONSTEXPR Number &
+ operator()(const TableIndices<rank_> &indices);
/**
* Return a @p const reference to the value referred to by the argument.
*/
- const Number &
- operator()(const TableIndices<rank_> &indices) const;
+ DEAL_II_CONSTEXPR const Number &
+ operator()(const TableIndices<rank_> &indices) const;
/**
* Access the elements of a row of this symmetric tensor. This function is
* called for constant tensors.
*/
- internal::SymmetricTensorAccessors::
+ constexpr internal::SymmetricTensorAccessors::
Accessor<rank_, dim, true, rank_ - 1, Number>
operator[](const unsigned int row) const;
* Access the elements of a row of this symmetric tensor. This function is
* called for non-constant tensors.
*/
- internal::SymmetricTensorAccessors::
+ DEAL_II_CONSTEXPR internal::SymmetricTensorAccessors::
Accessor<rank_, dim, false, rank_ - 1, Number>
operator[](const unsigned int row);
*
* Exactly the same as operator().
*/
- const Number &operator[](const TableIndices<rank_> &indices) const;
+ constexpr const Number &operator[](const TableIndices<rank_> &indices) const;
/**
* Return a read-write reference to the indicated element.
*
* Exactly the same as operator().
*/
- Number &operator[](const TableIndices<rank_> &indices);
+ DEAL_II_CONSTEXPR Number &operator[](const TableIndices<rank_> &indices);
/**
* Access to an element according to unrolled index. The function
* <tt>s.access_raw_entry(unrolled_index)</tt> does the same as
* <tt>s[s.unrolled_to_component_indices(i)]</tt>, but more efficiently.
*/
- const Number &
- access_raw_entry(const unsigned int unrolled_index) const;
+ DEAL_II_CONSTEXPR const Number &
+ access_raw_entry(const unsigned int unrolled_index) const;
/**
* Access to an element according to unrolled index. The function
* <tt>s.access_raw_entry(unrolled_index)</tt> does the same as
* <tt>s[s.unrolled_to_component_indices(i)]</tt>, but more efficiently.
*/
- Number &
- access_raw_entry(const unsigned int unrolled_index);
+ DEAL_II_CONSTEXPR Number &
+ access_raw_entry(const unsigned int unrolled_index);
/**
* Return the Frobenius-norm of a tensor, i.e. the square root of the sum of
* upper right as well as lower left entries, not just one of them, although
* they are equal for symmetric tensors).
*/
- typename numbers::NumberTraits<Number>::real_type
+ constexpr typename numbers::NumberTraits<Number>::real_type
norm() const;
/**
* <code>[0,n_independent_components)</code> the given entry in a symmetric
* tensor has.
*/
- static unsigned int
+ static constexpr unsigned int
component_to_unrolled_index(const TableIndices<rank_> &indices);
/**
* form of the tensor, return what set of indices $(k,l)$ (for rank-2
* tensors) or $(k,l,m,n)$ (for rank-4 tensors) corresponds to it.
*/
- static TableIndices<rank_>
+ static constexpr TableIndices<rank_>
unrolled_to_component_indices(const unsigned int i);
/**
* and indeed the state where all elements have a zero value is the state
* right after construction of such an object.
*/
- void
+ DEAL_II_CONSTEXPR void
clear();
/**
* Determine an estimate for the memory consumption (in bytes) of this
* object.
*/
- static std::size_t
+ static constexpr std::size_t
memory_consumption();
/**
* Make a few more functions friends.
*/
template <int dim2, typename Number2>
- friend Number2
- trace(const SymmetricTensor<2, dim2, Number2> &d);
+ friend DEAL_II_CONSTEXPR Number2
+ trace(const SymmetricTensor<2, dim2, Number2> &d);
template <int dim2, typename Number2>
- friend Number2
- determinant(const SymmetricTensor<2, dim2, Number2> &t);
+ friend DEAL_II_CONSTEXPR Number2
+ determinant(const SymmetricTensor<2, dim2, Number2> &t);
template <int dim2, typename Number2>
- friend SymmetricTensor<2, dim2, Number2>
- deviator(const SymmetricTensor<2, dim2, Number2> &t);
+ friend DEAL_II_CONSTEXPR SymmetricTensor<2, dim2, Number2>
+ deviator(const SymmetricTensor<2, dim2, Number2> &t);
template <int dim2, typename Number2>
- friend SymmetricTensor<2, dim2, Number2>
- unit_symmetric_tensor();
+ friend DEAL_II_CONSTEXPR SymmetricTensor<2, dim2, Number2>
+ unit_symmetric_tensor();
template <int dim2, typename Number2>
- friend SymmetricTensor<4, dim2, Number2>
- deviator_tensor();
+ friend DEAL_II_CONSTEXPR SymmetricTensor<4, dim2, Number2>
+ deviator_tensor();
template <int dim2, typename Number2>
- friend SymmetricTensor<4, dim2, Number2>
- identity_tensor();
+ friend DEAL_II_CONSTEXPR SymmetricTensor<4, dim2, Number2>
+ identity_tensor();
/**
namespace SymmetricTensorAccessors
{
template <int rank_, int dim, bool constness, int P, typename Number>
- Accessor<rank_, dim, constness, P, Number>::Accessor(
+ constexpr Accessor<rank_, dim, constness, P, Number>::Accessor(
tensor_type & tensor,
const TableIndices<rank_> &previous_indices)
: tensor(tensor)
template <int rank_, int dim, bool constness, int P, typename Number>
- Accessor<rank_, dim, constness, P - 1, Number>
+ DEAL_II_CONSTEXPR inline Accessor<rank_, dim, constness, P - 1, Number>
Accessor<rank_, dim, constness, P, Number>::
operator[](const unsigned int i)
{
template <int rank_, int dim, bool constness, int P, typename Number>
- Accessor<rank_, dim, constness, P - 1, Number>
+ constexpr Accessor<rank_, dim, constness, P - 1, Number>
Accessor<rank_, dim, constness, P, Number>::
operator[](const unsigned int i) const
{
template <int rank_, int dim, bool constness, typename Number>
- Accessor<rank_, dim, constness, 1, Number>::Accessor(
+ constexpr Accessor<rank_, dim, constness, 1, Number>::Accessor(
tensor_type & tensor,
const TableIndices<rank_> &previous_indices)
: tensor(tensor)
template <int rank_, int dim, bool constness, typename Number>
- typename Accessor<rank_, dim, constness, 1, Number>::reference
- Accessor<rank_, dim, constness, 1, Number>::
- operator[](const unsigned int i)
+ DEAL_II_CONSTEXPR inline
+ typename Accessor<rank_, dim, constness, 1, Number>::reference
+ Accessor<rank_, dim, constness, 1, Number>::
+ operator[](const unsigned int i)
{
return tensor(merge(previous_indices, i, rank_ - 1));
}
template <int rank_, int dim, bool constness, typename Number>
- typename Accessor<rank_, dim, constness, 1, Number>::reference
+ constexpr typename Accessor<rank_, dim, constness, 1, Number>::reference
Accessor<rank_, dim, constness, 1, Number>::
operator[](const unsigned int i) const
{
-template <int rank_, int dim, typename Number>
-inline SymmetricTensor<rank_, dim, Number>::SymmetricTensor()
-{
- // Some auto-differentiable numbers need explicit
- // zero initialization.
- for (unsigned int i = 0; i < base_tensor_type::dimension; ++i)
- data[i] = internal::NumberType<Number>::value(0.0);
-}
-
-
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline SymmetricTensor<rank_, dim, Number>::SymmetricTensor(
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim, Number>::SymmetricTensor(
const Tensor<2, dim, OtherNumber> &t)
{
Assert(rank == 2, ExcNotImplemented());
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline SymmetricTensor<rank_, dim, Number>::SymmetricTensor(
+constexpr SymmetricTensor<rank_, dim, Number>::SymmetricTensor(
const SymmetricTensor<rank_, dim, OtherNumber> &initializer)
-{
- for (unsigned int i = 0; i < base_tensor_type::dimension; ++i)
- data[i] =
- internal::NumberType<typename base_tensor_type::value_type>::value(
- initializer.data[i]);
-}
+ : data(initializer.data)
+{}
template <int rank_, int dim, typename Number>
-inline SymmetricTensor<rank_, dim, Number>::SymmetricTensor(
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim, Number>::SymmetricTensor(
const Number (&array)[n_independent_components])
: data(
*reinterpret_cast<const typename base_tensor_type::array_type *>(array))
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline SymmetricTensor<rank_, dim, Number> &
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim, Number> &
SymmetricTensor<rank_, dim, Number>::
operator=(const SymmetricTensor<rank_, dim, OtherNumber> &t)
{
- for (unsigned int i = 0; i < base_tensor_type::dimension; ++i)
- data[i] = t.data[i];
+ data = t.data;
return *this;
}
template <int rank_, int dim, typename Number>
-inline SymmetricTensor<rank_, dim, Number> &
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim, Number> &
SymmetricTensor<rank_, dim, Number>::operator=(const Number &d)
{
Assert(numbers::value_is_zero(d),
namespace SymmetricTensorImplementation
{
template <int dim, typename Number>
- inline DEAL_II_ALWAYS_INLINE dealii::Tensor<2, dim, Number>
- convert_to_tensor(const dealii::SymmetricTensor<2, dim, Number> &s)
+ DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ dealii::Tensor<2, dim, Number>
+ convert_to_tensor(const dealii::SymmetricTensor<2, dim, Number> &s)
{
dealii::Tensor<2, dim, Number> t;
template <int dim, typename Number>
- dealii::Tensor<4, dim, Number>
- convert_to_tensor(const dealii::SymmetricTensor<4, dim, Number> &st)
+ DEAL_II_CONSTEXPR dealii::Tensor<4, dim, Number>
+ convert_to_tensor(const dealii::SymmetricTensor<4, dim, Number> &st)
{
// utilize the symmetry properties of SymmetricTensor<4,dim>
// discussed in the class documentation to avoid accessing all
template <typename Number>
struct Inverse<2, 1, Number>
{
- static inline dealii::SymmetricTensor<2, 1, Number>
+ DEAL_II_CONSTEXPR static inline dealii::SymmetricTensor<2, 1, Number>
value(const dealii::SymmetricTensor<2, 1, Number> &t)
{
dealii::SymmetricTensor<2, 1, Number> tmp;
template <typename Number>
struct Inverse<2, 2, Number>
{
- static inline dealii::SymmetricTensor<2, 2, Number>
+ DEAL_II_CONSTEXPR static inline dealii::SymmetricTensor<2, 2, Number>
value(const dealii::SymmetricTensor<2, 2, Number> &t)
{
dealii::SymmetricTensor<2, 2, Number> tmp;
template <typename Number>
struct Inverse<2, 3, Number>
{
- static dealii::SymmetricTensor<2, 3, Number>
+ DEAL_II_CONSTEXPR static dealii::SymmetricTensor<2, 3, Number>
value(const dealii::SymmetricTensor<2, 3, Number> &t)
{
dealii::SymmetricTensor<2, 3, Number> tmp;
template <typename Number>
struct Inverse<4, 1, Number>
{
- static inline dealii::SymmetricTensor<4, 1, Number>
+ DEAL_II_CONSTEXPR static inline dealii::SymmetricTensor<4, 1, Number>
value(const dealii::SymmetricTensor<4, 1, Number> &t)
{
dealii::SymmetricTensor<4, 1, Number> tmp;
template <typename Number>
struct Inverse<4, 2, Number>
{
- static inline dealii::SymmetricTensor<4, 2, Number>
+ DEAL_II_CONSTEXPR static inline dealii::SymmetricTensor<4, 2, Number>
value(const dealii::SymmetricTensor<4, 2, Number> &t)
{
dealii::SymmetricTensor<4, 2, Number> tmp;
template <int rank_, int dim, typename Number>
-inline DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number>::
- operator Tensor<rank_, dim, Number>() const
+constexpr DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number>::
+ operator Tensor<rank_, dim, Number>() const
{
return internal::SymmetricTensorImplementation::convert_to_tensor(*this);
}
template <int rank_, int dim, typename Number>
-inline bool
+constexpr bool
SymmetricTensor<rank_, dim, Number>::
operator==(const SymmetricTensor<rank_, dim, Number> &t) const
{
template <int rank_, int dim, typename Number>
-inline bool
+constexpr bool
SymmetricTensor<rank_, dim, Number>::
operator!=(const SymmetricTensor<rank_, dim, Number> &t) const
{
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number> &
- SymmetricTensor<rank_, dim, Number>::
- operator+=(const SymmetricTensor<rank_, dim, OtherNumber> &t)
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ SymmetricTensor<rank_, dim, Number> &
+ SymmetricTensor<rank_, dim, Number>::
+ operator+=(const SymmetricTensor<rank_, dim, OtherNumber> &t)
{
data += t.data;
return *this;
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number> &
- SymmetricTensor<rank_, dim, Number>::
- operator-=(const SymmetricTensor<rank_, dim, OtherNumber> &t)
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ SymmetricTensor<rank_, dim, Number> &
+ SymmetricTensor<rank_, dim, Number>::
+ operator-=(const SymmetricTensor<rank_, dim, OtherNumber> &t)
{
data -= t.data;
return *this;
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number> &
-SymmetricTensor<rank_, dim, Number>::operator*=(const OtherNumber &d)
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ SymmetricTensor<rank_, dim, Number> &
+ SymmetricTensor<rank_, dim, Number>::operator*=(const OtherNumber &d)
{
data *= d;
return *this;
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number> &
-SymmetricTensor<rank_, dim, Number>::operator/=(const OtherNumber &d)
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ SymmetricTensor<rank_, dim, Number> &
+ SymmetricTensor<rank_, dim, Number>::operator/=(const OtherNumber &d)
{
data /= d;
return *this;
template <int rank_, int dim, typename Number>
-inline DEAL_II_ALWAYS_INLINE SymmetricTensor<rank_, dim, Number>
-SymmetricTensor<rank_, dim, Number>::operator-() const
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ SymmetricTensor<rank_, dim, Number>
+ SymmetricTensor<rank_, dim, Number>::operator-() const
{
SymmetricTensor tmp = *this;
tmp.data = -tmp.data;
template <int rank_, int dim, typename Number>
-inline DEAL_II_ALWAYS_INLINE void
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE void
SymmetricTensor<rank_, dim, Number>::clear()
{
data.clear();
template <int rank_, int dim, typename Number>
-inline std::size_t
+constexpr std::size_t
SymmetricTensor<rank_, dim, Number>::memory_consumption()
{
// all memory consists of statically allocated memory of the current
namespace internal
{
template <int dim, typename Number, typename OtherNumber = Number>
- inline DEAL_II_ALWAYS_INLINE typename SymmetricTensorAccessors::
- double_contraction_result<2, 2, dim, Number, OtherNumber>::type
- perform_double_contraction(
- const typename SymmetricTensorAccessors::StorageType<2, dim, Number>::
- base_tensor_type &data,
- const typename SymmetricTensorAccessors::
- StorageType<2, dim, OtherNumber>::base_tensor_type &sdata)
+ DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ typename SymmetricTensorAccessors::
+ double_contraction_result<2, 2, dim, Number, OtherNumber>::type
+ perform_double_contraction(
+ const typename SymmetricTensorAccessors::StorageType<2, dim, Number>::
+ base_tensor_type &data,
+ const typename SymmetricTensorAccessors::
+ StorageType<2, dim, OtherNumber>::base_tensor_type &sdata)
{
using result_type = typename SymmetricTensorAccessors::
double_contraction_result<2, 2, dim, Number, OtherNumber>::type;
template <int dim, typename Number, typename OtherNumber = Number>
- inline typename SymmetricTensorAccessors::
+ DEAL_II_CONSTEXPR inline typename SymmetricTensorAccessors::
double_contraction_result<4, 2, dim, Number, OtherNumber>::type
perform_double_contraction(
const typename SymmetricTensorAccessors::StorageType<4, dim, Number>::
template <int dim, typename Number, typename OtherNumber = Number>
- inline typename SymmetricTensorAccessors::StorageType<
+ DEAL_II_CONSTEXPR inline typename SymmetricTensorAccessors::StorageType<
2,
dim,
typename SymmetricTensorAccessors::
template <int dim, typename Number, typename OtherNumber = Number>
- inline typename SymmetricTensorAccessors::StorageType<
+ DEAL_II_CONSTEXPR inline typename SymmetricTensorAccessors::StorageType<
4,
dim,
typename SymmetricTensorAccessors::
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE typename internal::SymmetricTensorAccessors::
- double_contraction_result<rank_, 2, dim, Number, OtherNumber>::type
- SymmetricTensor<rank_, dim, Number>::
- operator*(const SymmetricTensor<2, dim, OtherNumber> &s) const
+DEAL_II_CONSTEXPR inline DEAL_II_ALWAYS_INLINE
+ typename internal::SymmetricTensorAccessors::
+ double_contraction_result<rank_, 2, dim, Number, OtherNumber>::type
+ SymmetricTensor<rank_, dim, Number>::
+ operator*(const SymmetricTensor<2, dim, OtherNumber> &s) const
{
// need to have two different function calls
// because a scalar and rank-2 tensor are not
template <int rank_, int dim, typename Number>
template <typename OtherNumber>
-inline typename internal::SymmetricTensorAccessors::
+DEAL_II_CONSTEXPR inline typename internal::SymmetricTensorAccessors::
double_contraction_result<rank_, 4, dim, Number, OtherNumber>::type
SymmetricTensor<rank_, dim, Number>::
operator*(const SymmetricTensor<4, dim, OtherNumber> &s) const
// into a separate namespace
namespace internal
{
+ // The variables within this struct will be referenced in the next functions.
+ // It is a workaround that allows returning a reference to a static variable
+ // while allowing constexpr evaluation of the function.
+ // It has to be defined outside the function because constexpr functions
+ // cannot define static variables.
+ // A similar struct has also been defined in tensor.h
+ template <typename Type>
+ struct Uninitialized
+ {
+ static Type value;
+ };
+
+ template <typename Type>
+ Type Uninitialized<Type>::value;
+
template <int dim, typename Number>
- inline Number &
+ DEAL_II_CONSTEXPR inline Number &
symmetric_tensor_access(const TableIndices<2> &indices,
typename SymmetricTensorAccessors::
StorageType<2, dim, Number>::base_tensor_type &data)
}
}
- static Number dummy_but_referenceable = Number();
- return dummy_but_referenceable;
+ // The code should never reach there.
+ // Returns a dummy reference to a dummy variable just to make the
+ // compiler happy.
+ return Uninitialized<Number>::value;
}
template <int dim, typename Number>
- inline const Number &
+ DEAL_II_CONSTEXPR inline const Number &
symmetric_tensor_access(const TableIndices<2> &indices,
const typename SymmetricTensorAccessors::
StorageType<2, dim, Number>::base_tensor_type &data)
default:
// to do the rest, sort our indices before comparing
{
- TableIndices<2> sorted_indices(indices);
- sorted_indices.sort();
-
+ TableIndices<2> sorted_indices(std::min(indices[0], indices[1]),
+ std::max(indices[0], indices[1]));
for (unsigned int d = 0, c = 0; d < dim; ++d)
for (unsigned int e = d + 1; e < dim; ++e, ++c)
if ((sorted_indices[0] == d) && (sorted_indices[1] == e))
}
}
- static Number dummy_but_referenceable = Number();
- return dummy_but_referenceable;
+ // The code should never reach there.
+ // Returns a dummy reference to a dummy variable just to make the
+ // compiler happy.
+ return Uninitialized<Number>::value;
}
template <int dim, typename Number>
- inline Number &
+ DEAL_II_CONSTEXPR inline Number &
symmetric_tensor_access(const TableIndices<4> &indices,
typename SymmetricTensorAccessors::
StorageType<4, dim, Number>::base_tensor_type &data)
return data[0][0];
case 2:
- // each entry of the tensor can be
- // thought of as an entry in a
- // matrix that maps the rolled-out
- // rank-2 tensors into rolled-out
- // rank-2 tensors. this is the
- // format in which we store rank-4
- // tensors. determine which
- // position the present entry is
+ // each entry of the tensor can be thought of as an entry in a
+ // matrix that maps the rolled-out rank-2 tensors into rolled-out
+ // rank-2 tensors. this is the format in which we store rank-4
+ // tensors. determine which position the present entry is
// stored in
{
- unsigned int base_index[2];
- if ((indices[0] == 0) && (indices[1] == 0))
- base_index[0] = 0;
- else if ((indices[0] == 1) && (indices[1] == 1))
- base_index[0] = 1;
- else
- base_index[0] = 2;
-
- if ((indices[2] == 0) && (indices[3] == 0))
- base_index[1] = 0;
- else if ((indices[2] == 1) && (indices[3] == 1))
- base_index[1] = 1;
- else
- base_index[1] = 2;
-
- return data[base_index[0]][base_index[1]];
+ constexpr std::size_t base_index[2][2] = {{0, 2}, {2, 1}};
+ return data[base_index[indices[0]][indices[1]]]
+ [base_index[indices[2]][indices[3]]];
}
-
case 3:
- // each entry of the tensor can be
- // thought of as an entry in a
- // matrix that maps the rolled-out
- // rank-2 tensors into rolled-out
- // rank-2 tensors. this is the
- // format in which we store rank-4
- // tensors. determine which
- // position the present entry is
+ // each entry of the tensor can be thought of as an entry in a
+ // matrix that maps the rolled-out rank-2 tensors into rolled-out
+ // rank-2 tensors. this is the format in which we store rank-4
+ // tensors. determine which position the present entry is
// stored in
{
- unsigned int base_index[2];
- if ((indices[0] == 0) && (indices[1] == 0))
- base_index[0] = 0;
- else if ((indices[0] == 1) && (indices[1] == 1))
- base_index[0] = 1;
- else if ((indices[0] == 2) && (indices[1] == 2))
- base_index[0] = 2;
- else if (((indices[0] == 0) && (indices[1] == 1)) ||
- ((indices[0] == 1) && (indices[1] == 0)))
- base_index[0] = 3;
- else if (((indices[0] == 0) && (indices[1] == 2)) ||
- ((indices[0] == 2) && (indices[1] == 0)))
- base_index[0] = 4;
- else
- {
- Assert(((indices[0] == 1) && (indices[1] == 2)) ||
- ((indices[0] == 2) && (indices[1] == 1)),
- ExcInternalError());
- base_index[0] = 5;
- }
-
- if ((indices[2] == 0) && (indices[3] == 0))
- base_index[1] = 0;
- else if ((indices[2] == 1) && (indices[3] == 1))
- base_index[1] = 1;
- else if ((indices[2] == 2) && (indices[3] == 2))
- base_index[1] = 2;
- else if (((indices[2] == 0) && (indices[3] == 1)) ||
- ((indices[2] == 1) && (indices[3] == 0)))
- base_index[1] = 3;
- else if (((indices[2] == 0) && (indices[3] == 2)) ||
- ((indices[2] == 2) && (indices[3] == 0)))
- base_index[1] = 4;
- else
- {
- Assert(((indices[2] == 1) && (indices[3] == 2)) ||
- ((indices[2] == 2) && (indices[3] == 1)),
- ExcInternalError());
- base_index[1] = 5;
- }
-
- return data[base_index[0]][base_index[1]];
+ constexpr std::size_t base_index[3][3] = {{0, 3, 4},
+ {3, 1, 5},
+ {4, 5, 2}};
+ return data[base_index[indices[0]][indices[1]]]
+ [base_index[indices[2]][indices[3]]];
}
default:
Assert(false, ExcNotImplemented());
}
- static Number dummy;
- return dummy;
+ // The code should never reach there.
+ // Returns a dummy reference to a dummy variable just to make the
+ // compiler happy.
+ return Uninitialized<Number>::value;
}
template <int dim, typename Number>
- inline const Number &
+ DEAL_II_CONSTEXPR inline const Number &
symmetric_tensor_access(const TableIndices<4> &indices,
const typename SymmetricTensorAccessors::
StorageType<4, dim, Number>::base_tensor_type &data)
return data[0][0];
case 2:
- // each entry of the tensor can be
- // thought of as an entry in a
- // matrix that maps the rolled-out
- // rank-2 tensors into rolled-out
- // rank-2 tensors. this is the
- // format in which we store rank-4
- // tensors. determine which
- // position the present entry is
+ // each entry of the tensor can be thought of as an entry in a
+ // matrix that maps the rolled-out rank-2 tensors into rolled-out
+ // rank-2 tensors. this is the format in which we store rank-4
+ // tensors. determine which position the present entry is
// stored in
{
- unsigned int base_index[2];
- if ((indices[0] == 0) && (indices[1] == 0))
- base_index[0] = 0;
- else if ((indices[0] == 1) && (indices[1] == 1))
- base_index[0] = 1;
- else
- base_index[0] = 2;
-
- if ((indices[2] == 0) && (indices[3] == 0))
- base_index[1] = 0;
- else if ((indices[2] == 1) && (indices[3] == 1))
- base_index[1] = 1;
- else
- base_index[1] = 2;
-
- return data[base_index[0]][base_index[1]];
+ constexpr std::size_t base_index[2][2] = {{0, 2}, {2, 1}};
+ return data[base_index[indices[0]][indices[1]]]
+ [base_index[indices[2]][indices[3]]];
}
-
case 3:
- // each entry of the tensor can be
- // thought of as an entry in a
- // matrix that maps the rolled-out
- // rank-2 tensors into rolled-out
- // rank-2 tensors. this is the
- // format in which we store rank-4
- // tensors. determine which
- // position the present entry is
+ // each entry of the tensor can be thought of as an entry in a
+ // matrix that maps the rolled-out rank-2 tensors into rolled-out
+ // rank-2 tensors. this is the format in which we store rank-4
+ // tensors. determine which position the present entry is
// stored in
{
- unsigned int base_index[2];
- if ((indices[0] == 0) && (indices[1] == 0))
- base_index[0] = 0;
- else if ((indices[0] == 1) && (indices[1] == 1))
- base_index[0] = 1;
- else if ((indices[0] == 2) && (indices[1] == 2))
- base_index[0] = 2;
- else if (((indices[0] == 0) && (indices[1] == 1)) ||
- ((indices[0] == 1) && (indices[1] == 0)))
- base_index[0] = 3;
- else if (((indices[0] == 0) && (indices[1] == 2)) ||
- ((indices[0] == 2) && (indices[1] == 0)))
- base_index[0] = 4;
- else
- {
- Assert(((indices[0] == 1) && (indices[1] == 2)) ||
- ((indices[0] == 2) && (indices[1] == 1)),
- ExcInternalError());
- base_index[0] = 5;
- }
-
- if ((indices[2] == 0) && (indices[3] == 0))
- base_index[1] = 0;
- else if ((indices[2] == 1) && (indices[3] == 1))
- base_index[1] = 1;
- else if ((indices[2] == 2) && (indices[3] == 2))
- base_index[1] = 2;
- else if (((indices[2] == 0) && (indices[3] == 1)) ||
- ((indices[2] == 1) && (indices[3] == 0)))
- base_index[1] = 3;
- else if (((indices[2] == 0) && (indices[3] == 2)) ||
- ((indices[2] == 2) && (indices[3] == 0)))
- base_index[1] = 4;
- else
- {
- Assert(((indices[2] == 1) && (indices[3] == 2)) ||
- ((indices[2] == 2) && (indices[3] == 1)),
- ExcInternalError());
- base_index[1] = 5;
- }
-
- return data[base_index[0]][base_index[1]];
+ constexpr std::size_t base_index[3][3] = {{0, 3, 4},
+ {3, 1, 5},
+ {4, 5, 2}};
+ return data[base_index[indices[0]][indices[1]]]
+ [base_index[indices[2]][indices[3]]];
}
default:
Assert(false, ExcNotImplemented());
}
- static Number dummy;
- return dummy;
+ // The code should never reach there.
+ // Returns a dummy reference to a dummy variable just to make the
+ // compiler happy.
+ return Uninitialized<Number>::value;
}
} // end of namespace internal
template <int rank_, int dim, typename Number>
-inline Number &
+DEAL_II_CONSTEXPR inline Number &
SymmetricTensor<rank_, dim, Number>::
operator()(const TableIndices<rank_> &indices)
{
template <int rank_, int dim, typename Number>
-inline const Number &
+DEAL_II_CONSTEXPR inline const Number &
SymmetricTensor<rank_, dim, Number>::
operator()(const TableIndices<rank_> &indices) const
{
namespace SymmetricTensorImplementation
{
template <int rank_>
- TableIndices<rank_>
+ constexpr TableIndices<rank_>
get_partially_filled_indices(const unsigned int row,
const std::integral_constant<int, 2> &)
{
template <int rank_>
- TableIndices<rank_>
+ constexpr TableIndices<rank_>
get_partially_filled_indices(const unsigned int row,
const std::integral_constant<int, 4> &)
{
template <int rank_, int dim, typename Number>
-internal::SymmetricTensorAccessors::
+constexpr internal::SymmetricTensorAccessors::
Accessor<rank_, dim, true, rank_ - 1, Number>
SymmetricTensor<rank_, dim, Number>::
operator[](const unsigned int row) const
template <int rank_, int dim, typename Number>
-internal::SymmetricTensorAccessors::
+DEAL_II_CONSTEXPR inline internal::SymmetricTensorAccessors::
Accessor<rank_, dim, false, rank_ - 1, Number>
SymmetricTensor<rank_, dim, Number>::operator[](const unsigned int row)
{
template <int rank_, int dim, typename Number>
-inline const Number &SymmetricTensor<rank_, dim, Number>::
- operator[](const TableIndices<rank_> &indices) const
+constexpr const Number &SymmetricTensor<rank_, dim, Number>::
+ operator[](const TableIndices<rank_> &indices) const
{
return operator()(indices);
}
template <int rank_, int dim, typename Number>
-inline Number &SymmetricTensor<rank_, dim, Number>::
- operator[](const TableIndices<rank_> &indices)
+DEAL_II_CONSTEXPR inline Number &SymmetricTensor<rank_, dim, Number>::
+ operator[](const TableIndices<rank_> &indices)
{
return operator()(indices);
}
template <int rank_, int dim, typename Number>
-inline Number *
+DEAL_II_CONSTEXPR inline Number *
SymmetricTensor<rank_, dim, Number>::begin_raw()
{
return std::addressof(this->access_raw_entry(0));
template <int rank_, int dim, typename Number>
-inline const Number *
+constexpr const Number *
SymmetricTensor<rank_, dim, Number>::begin_raw() const
{
return std::addressof(this->access_raw_entry(0));
template <int rank_, int dim, typename Number>
-inline Number *
+DEAL_II_CONSTEXPR inline Number *
SymmetricTensor<rank_, dim, Number>::end_raw()
{
return begin_raw() + n_independent_components;
template <int rank_, int dim, typename Number>
-inline const Number *
+constexpr const Number *
SymmetricTensor<rank_, dim, Number>::end_raw() const
{
return begin_raw() + n_independent_components;
namespace SymmetricTensorImplementation
{
template <int dim, typename Number>
- unsigned int
+ constexpr unsigned int
entry_to_indices(const dealii::SymmetricTensor<2, dim, Number> &,
const unsigned int index)
{
template <int dim, typename Number>
- dealii::TableIndices<2>
+ constexpr dealii::TableIndices<2>
entry_to_indices(const dealii::SymmetricTensor<4, dim, Number> &,
const unsigned int index)
{
template <int rank_, int dim, typename Number>
-inline const Number &
+DEAL_II_CONSTEXPR inline const Number &
SymmetricTensor<rank_, dim, Number>::access_raw_entry(
const unsigned int index) const
{
template <int rank_, int dim, typename Number>
-inline Number &
+DEAL_II_CONSTEXPR inline Number &
SymmetricTensor<rank_, dim, Number>::access_raw_entry(const unsigned int index)
{
AssertIndexRange(index, n_independent_components);
namespace internal
{
template <int dim, typename Number>
- inline typename numbers::NumberTraits<Number>::real_type
+ DEAL_II_CONSTEXPR inline typename numbers::NumberTraits<Number>::real_type
compute_norm(const typename SymmetricTensorAccessors::
StorageType<2, dim, Number>::base_tensor_type &data)
{
template <int dim, typename Number>
- inline typename numbers::NumberTraits<Number>::real_type
+ DEAL_II_CONSTEXPR inline typename numbers::NumberTraits<Number>::real_type
compute_norm(const typename SymmetricTensorAccessors::
StorageType<4, dim, Number>::base_tensor_type &data)
{
template <int rank_, int dim, typename Number>
-inline typename numbers::NumberTraits<Number>::real_type
+constexpr typename numbers::NumberTraits<Number>::real_type
SymmetricTensor<rank_, dim, Number>::norm() const
{
return internal::compute_norm<dim, Number>(data);
//
// this function is for rank-2 tensors
template <int dim>
- inline unsigned int
+ DEAL_II_CONSTEXPR inline unsigned int
component_to_unrolled_index(const TableIndices<2> &indices)
{
Assert(indices[0] < dim, ExcIndexRange(indices[0], 0, dim));
case 2:
{
- static const unsigned int table[2][2] = {{0, 2}, {2, 1}};
+ constexpr unsigned int table[2][2] = {{0, 2}, {2, 1}};
return table[indices[0]][indices[1]];
}
case 3:
{
- static const unsigned int table[3][3] = {{0, 3, 4},
- {3, 1, 5},
- {4, 5, 2}};
+ constexpr unsigned int table[3][3] = {{0, 3, 4},
+ {3, 1, 5},
+ {4, 5, 2}};
return table[indices[0]][indices[1]];
}
case 4:
{
- static const unsigned int table[4][4] = {{0, 4, 5, 6},
- {4, 1, 7, 8},
- {5, 7, 2, 9},
- {6, 8, 9, 3}};
+ constexpr unsigned int table[4][4] = {{0, 4, 5, 6},
+ {4, 1, 7, 8},
+ {5, 7, 2, 9},
+ {6, 8, 9, 3}};
return table[indices[0]][indices[1]];
}
// this function is for tensors of ranks not already handled
// above
template <int dim, int rank_>
- inline unsigned int
+ DEAL_II_CONSTEXPR inline unsigned int
component_to_unrolled_index(const TableIndices<rank_> &indices)
{
(void)indices;
template <int rank_, int dim, typename Number>
-inline unsigned int
+constexpr unsigned int
SymmetricTensor<rank_, dim, Number>::component_to_unrolled_index(
const TableIndices<rank_> &indices)
{
//
// this function is for rank-2 tensors
template <int dim>
- inline TableIndices<2>
+ DEAL_II_CONSTEXPR inline TableIndices<2>
unrolled_to_component_indices(const unsigned int i,
const std::integral_constant<int, 2> &)
{
// this function is for tensors of a rank not already handled
// above
template <int dim, int rank_>
- inline TableIndices<rank_>
+ DEAL_II_CONSTEXPR inline TableIndices<rank_>
unrolled_to_component_indices(const unsigned int i,
const std::integral_constant<int, rank_> &)
{
} // namespace internal
template <int rank_, int dim, typename Number>
-inline TableIndices<rank_>
+constexpr TableIndices<rank_>
SymmetricTensor<rank_, dim, Number>::unrolled_to_component_indices(
const unsigned int i)
{
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline SymmetricTensor<rank_,
- dim,
- typename ProductType<Number, OtherNumber>::type>
+DEAL_II_CONSTEXPR inline SymmetricTensor<
+ rank_,
+ dim,
+ typename ProductType<Number, OtherNumber>::type>
operator+(const SymmetricTensor<rank_, dim, Number> & left,
const SymmetricTensor<rank_, dim, OtherNumber> &right)
{
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline SymmetricTensor<rank_,
- dim,
- typename ProductType<Number, OtherNumber>::type>
+DEAL_II_CONSTEXPR inline SymmetricTensor<
+ rank_,
+ dim,
+ typename ProductType<Number, OtherNumber>::type>
operator-(const SymmetricTensor<rank_, dim, Number> & left,
const SymmetricTensor<rank_, dim, OtherNumber> &right)
{
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
+constexpr Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
operator+(const SymmetricTensor<rank_, dim, Number> &left,
const Tensor<rank_, dim, OtherNumber> & right)
{
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
+constexpr Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
operator+(const Tensor<rank_, dim, Number> & left,
const SymmetricTensor<rank_, dim, OtherNumber> &right)
{
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
+constexpr Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
operator-(const SymmetricTensor<rank_, dim, Number> &left,
const Tensor<rank_, dim, OtherNumber> & right)
{
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
+constexpr Tensor<rank_, dim, typename ProductType<Number, OtherNumber>::type>
operator-(const Tensor<rank_, dim, Number> & left,
const SymmetricTensor<rank_, dim, OtherNumber> &right)
{
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline Number
+DEAL_II_CONSTEXPR inline Number
determinant(const SymmetricTensor<2, dim, Number> &t)
{
switch (dim)
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline Number
+constexpr Number
third_invariant(const SymmetricTensor<2, dim, Number> &t)
{
return determinant(t);
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-Number
-trace(const SymmetricTensor<2, dim, Number> &d)
+DEAL_II_CONSTEXPR Number
+ trace(const SymmetricTensor<2, dim, Number> &d)
{
Number t = d.data[0];
for (unsigned int i = 1; i < dim; ++i)
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline Number
+constexpr Number
first_invariant(const SymmetricTensor<2, dim, Number> &t)
{
return trace(t);
* @author Wolfgang Bangerth, 2005, 2010
*/
template <typename Number>
-inline Number
+constexpr Number
second_invariant(const SymmetricTensor<2, 1, Number> &)
{
return internal::NumberType<Number>::value(0.0);
* @author Wolfgang Bangerth, 2005, 2010
*/
template <typename Number>
-inline Number
+constexpr Number
second_invariant(const SymmetricTensor<2, 2, Number> &t)
{
return t[0][0] * t[1][1] - t[0][1] * t[0][1];
* @author Wolfgang Bangerth, 2005, 2010
*/
template <typename Number>
-inline Number
+constexpr Number
second_invariant(const SymmetricTensor<2, 3, Number> &t)
{
return (t[0][0] * t[1][1] + t[1][1] * t[2][2] + t[2][2] * t[0][0] -
* @author Wolfgang Bangerth, 2005
*/
template <int rank_, int dim, typename Number>
-inline SymmetricTensor<rank_, dim, Number>
+constexpr SymmetricTensor<rank_, dim, Number>
transpose(const SymmetricTensor<rank_, dim, Number> &t)
{
return t;
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<2, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<2, dim, Number>
deviator(const SymmetricTensor<2, dim, Number> &t)
{
SymmetricTensor<2, dim, Number> tmp = t;
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<2, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<2, dim, Number>
unit_symmetric_tensor()
{
// create a default constructed matrix filled with
* @author Wolfgang Bangerth, 2005
*/
template <int dim>
-inline SymmetricTensor<2, dim>
+DEAL_II_CONSTEXPR inline SymmetricTensor<2, dim>
unit_symmetric_tensor()
{
return unit_symmetric_tensor<dim, double>();
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<4, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<4, dim, Number>
deviator_tensor()
{
SymmetricTensor<4, dim, Number> tmp;
* @author Wolfgang Bangerth, 2005
*/
template <int dim>
-inline SymmetricTensor<4, dim>
+DEAL_II_CONSTEXPR inline SymmetricTensor<4, dim>
deviator_tensor()
{
return deviator_tensor<dim, double>();
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<4, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<4, dim, Number>
identity_tensor()
{
SymmetricTensor<4, dim, Number> tmp;
* @author Wolfgang Bangerth, 2005
*/
template <int dim>
-inline SymmetricTensor<4, dim>
+DEAL_II_CONSTEXPR inline SymmetricTensor<4, dim>
identity_tensor()
{
return identity_tensor<dim, double>();
* @author Jean-Paul Pelteret, 2016
*/
template <int dim, typename Number>
-inline SymmetricTensor<2, dim, Number>
+constexpr SymmetricTensor<2, dim, Number>
invert(const SymmetricTensor<2, dim, Number> &t)
{
return internal::SymmetricTensorImplementation::Inverse<2, dim, Number>::
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<4, dim, Number>
+constexpr SymmetricTensor<4, dim, Number>
invert(const SymmetricTensor<4, dim, Number> &t)
{
return internal::SymmetricTensorImplementation::Inverse<4, dim, Number>::
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<4, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<4, dim, Number>
outer_product(const SymmetricTensor<2, dim, Number> &t1,
const SymmetricTensor<2, dim, Number> &t2)
{
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number>
-inline SymmetricTensor<2, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<2, dim, Number>
symmetrize(const Tensor<2, dim, Number> &t)
{
Number array[(dim * dim + dim) / 2];
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number>
-inline SymmetricTensor<rank_, dim, Number>
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim, Number>
operator*(const SymmetricTensor<rank_, dim, Number> &t, const Number &factor)
{
SymmetricTensor<rank_, dim, Number> tt = t;
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number>
-inline SymmetricTensor<rank_, dim, Number>
+constexpr SymmetricTensor<rank_, dim, Number>
operator*(const Number &factor, const SymmetricTensor<rank_, dim, Number> &t)
{
// simply forward to the other operator
* @relatesalso EnableIfScalar
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline SymmetricTensor<
+DEAL_II_CONSTEXPR inline SymmetricTensor<
rank_,
dim,
typename ProductType<Number,
* @relatesalso EnableIfScalar
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline SymmetricTensor<
+DEAL_II_CONSTEXPR inline SymmetricTensor<
rank_,
dim,
typename ProductType<OtherNumber,
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim, typename Number, typename OtherNumber>
-inline SymmetricTensor<
+DEAL_II_CONSTEXPR inline SymmetricTensor<
rank_,
dim,
typename ProductType<Number,
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim>
-inline SymmetricTensor<rank_, dim>
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim>
operator*(const SymmetricTensor<rank_, dim> &t, const double factor)
{
SymmetricTensor<rank_, dim> tt = t;
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim>
-inline SymmetricTensor<rank_, dim>
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim>
operator*(const double factor, const SymmetricTensor<rank_, dim> &t)
{
SymmetricTensor<rank_, dim> tt = t;
* @relatesalso SymmetricTensor
*/
template <int rank_, int dim>
-inline SymmetricTensor<rank_, dim>
+DEAL_II_CONSTEXPR inline SymmetricTensor<rank_, dim>
operator/(const SymmetricTensor<rank_, dim> &t, const double factor)
{
SymmetricTensor<rank_, dim> tt = t;
* @relatesalso SymmetricTensor
*/
template <int dim, typename Number, typename OtherNumber>
-inline typename ProductType<Number, OtherNumber>::type
+constexpr typename ProductType<Number, OtherNumber>::type
scalar_product(const SymmetricTensor<2, dim, Number> & t1,
const SymmetricTensor<2, dim, OtherNumber> &t2)
{
* @relatesalso Tensor @relatesalso SymmetricTensor
*/
template <int dim, typename Number, typename OtherNumber>
-inline typename ProductType<Number, OtherNumber>::type
+DEAL_II_CONSTEXPR inline typename ProductType<Number, OtherNumber>::type
scalar_product(const SymmetricTensor<2, dim, Number> &t1,
const Tensor<2, dim, OtherNumber> & t2)
{
* @relatesalso Tensor @relatesalso SymmetricTensor
*/
template <int dim, typename Number, typename OtherNumber>
-inline typename ProductType<Number, OtherNumber>::type
+constexpr typename ProductType<Number, OtherNumber>::type
scalar_product(const Tensor<2, dim, Number> & t1,
const SymmetricTensor<2, dim, OtherNumber> &t2)
{
* @author Wolfgang Bangerth, 2005
*/
template <typename Number, typename OtherNumber>
-inline void double_contract(
+DEAL_II_CONSTEXPR inline void double_contract(
SymmetricTensor<2, 1, typename ProductType<Number, OtherNumber>::type> &tmp,
const SymmetricTensor<4, 1, Number> & t,
const SymmetricTensor<2, 1, OtherNumber> & s)
* @author Wolfgang Bangerth, 2005
*/
template <typename Number, typename OtherNumber>
-inline void double_contract(
+DEAL_II_CONSTEXPR inline void double_contract(
SymmetricTensor<2, 1, typename ProductType<Number, OtherNumber>::type> &tmp,
const SymmetricTensor<2, 1, Number> & s,
const SymmetricTensor<4, 1, OtherNumber> & t)
* @relatesalso SymmetricTensor @author Wolfgang Bangerth, 2005
*/
template <typename Number, typename OtherNumber>
-inline void double_contract(
+DEAL_II_CONSTEXPR inline void double_contract(
SymmetricTensor<2, 2, typename ProductType<Number, OtherNumber>::type> &tmp,
const SymmetricTensor<4, 2, Number> & t,
const SymmetricTensor<2, 2, OtherNumber> & s)
* @author Wolfgang Bangerth, 2005
*/
template <typename Number, typename OtherNumber>
-inline void double_contract(
+DEAL_II_CONSTEXPR inline void double_contract(
SymmetricTensor<2, 2, typename ProductType<Number, OtherNumber>::type> &tmp,
const SymmetricTensor<2, 2, Number> & s,
const SymmetricTensor<4, 2, OtherNumber> & t)
* @author Wolfgang Bangerth, 2005
*/
template <typename Number, typename OtherNumber>
-inline void double_contract(
+DEAL_II_CONSTEXPR inline void double_contract(
SymmetricTensor<2, 3, typename ProductType<Number, OtherNumber>::type> &tmp,
const SymmetricTensor<4, 3, Number> & t,
const SymmetricTensor<2, 3, OtherNumber> & s)
* @author Wolfgang Bangerth, 2005
*/
template <typename Number, typename OtherNumber>
-inline void double_contract(
+DEAL_II_CONSTEXPR inline void double_contract(
SymmetricTensor<2, 3, typename ProductType<Number, OtherNumber>::type> &tmp,
const SymmetricTensor<2, 3, Number> & s,
const SymmetricTensor<4, 3, OtherNumber> & t)
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number, typename OtherNumber>
-Tensor<1, dim, typename ProductType<Number, OtherNumber>::type>
-operator*(const SymmetricTensor<2, dim, Number> &src1,
- const Tensor<1, dim, OtherNumber> & src2)
+DEAL_II_CONSTEXPR
+ Tensor<1, dim, typename ProductType<Number, OtherNumber>::type>
+ operator*(const SymmetricTensor<2, dim, Number> &src1,
+ const Tensor<1, dim, OtherNumber> & src2)
{
Tensor<1, dim, typename ProductType<Number, OtherNumber>::type> dest;
for (unsigned int i = 0; i < dim; ++i)
* @author Wolfgang Bangerth, 2005
*/
template <int dim, typename Number, typename OtherNumber>
-Tensor<1, dim, typename ProductType<Number, OtherNumber>::type>
+constexpr Tensor<1, dim, typename ProductType<Number, OtherNumber>::type>
operator*(const Tensor<1, dim, Number> & src1,
const SymmetricTensor<2, dim, OtherNumber> &src2)
{
int dim,
typename Number,
typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE
+constexpr DEAL_II_ALWAYS_INLINE
typename Tensor<rank_1 + rank_2 - 2,
dim,
typename ProductType<Number, OtherNumber>::type>::tensor_type
operator*(const Tensor<rank_1, dim, Number> & src1,
- const SymmetricTensor<rank_2, dim, OtherNumber> &src2s)
+ const SymmetricTensor<rank_2, dim, OtherNumber> &src2)
{
- const Tensor<rank_2, dim, OtherNumber> src2(src2s);
- return src1 * src2;
+ return src1 * Tensor<rank_2, dim, OtherNumber>(src2);
}
int dim,
typename Number,
typename OtherNumber>
-inline DEAL_II_ALWAYS_INLINE
+constexpr DEAL_II_ALWAYS_INLINE
typename Tensor<rank_1 + rank_2 - 2,
dim,
typename ProductType<Number, OtherNumber>::type>::tensor_type
- operator*(const SymmetricTensor<rank_1, dim, Number> &src1s,
+ operator*(const SymmetricTensor<rank_1, dim, Number> &src1,
const Tensor<rank_2, dim, OtherNumber> & src2)
{
- const Tensor<rank_2, dim, OtherNumber> src1(src1s);
- return src1 * src2;
+ return Tensor<rank_2, dim, OtherNumber>(src1) * src2;
}