VoigtArray.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 Equinor ASA
3
4 This file is part of the Open Porous Media Project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_UTIL_VOIGT_ARRAY_HPP
21#define OPM_UTIL_VOIGT_ARRAY_HPP
22
23#include <algorithm>
24#include <array>
25#include <cstddef>
26#include <initializer_list>
27#include <type_traits>
28#include <vector>
29
30namespace Opm {
31
32enum class VoigtIndex {
33 XX = 0, XY = 5, XZ = 4,
34 YX = XY, YY = 1, YZ = 3,
35 ZX = XZ, ZY = YZ, ZZ = 2,
36};
37
38template<class T>
40{
41public:
42 static constexpr auto indices = std::array<VoigtIndex, 9>{ // NVCC needs templates specified
52 };
53
54 static constexpr auto unique_indices = std::array<VoigtIndex, 6>{ // NVCC needs templates specified
61 };
62
63 static constexpr auto diag_indices = std::array<VoigtIndex, 3>{ // NVCC needs templates specified
67 };
68
69 VoigtContainer() = default;
70
71 template<class Array>
72 VoigtContainer(const Array& array);
73
74 VoigtContainer(std::initializer_list<T> value)
75 {
76 std::copy_n(value.begin(),
77 std::min(data_.size(), value.size()),
78 data_.begin());
79 }
80
81 const T& operator[](const VoigtIndex idx) const
82 { return data_[static_cast<std::underlying_type_t<VoigtIndex>>(idx)]; }
83
84 T& operator [](const VoigtIndex idx)
85 { return data_[static_cast<std::underlying_type_t<VoigtIndex>>(idx)]; }
86
87 constexpr std::size_t size() const { return data_.size(); }
88
89protected:
90 std::array<T, 6> data_{};
91};
92
93template<class Scalar>
94class VoigtArray : public VoigtContainer<std::vector<Scalar>>
95{
96public:
97 VoigtArray() = default;
98 explicit VoigtArray(const std::size_t size);
99
100 void resize(const std::size_t size);
101
102 Scalar operator()(const VoigtIndex idx, const std::size_t i) const;
103 Scalar& operator()(const VoigtIndex idx, const std::size_t i);
104
105 void assign(const std::size_t i, const VoigtContainer<Scalar>& array);
106};
107
108} // namespace Opm
109
110#endif // OPM_UTIL_VOIGT_ARRAY_HPP
Definition: VoigtArray.hpp:95
Scalar & operator()(const VoigtIndex idx, const std::size_t i)
VoigtArray()=default
VoigtArray(const std::size_t size)
void assign(const std::size_t i, const VoigtContainer< Scalar > &array)
Scalar operator()(const VoigtIndex idx, const std::size_t i) const
void resize(const std::size_t size)
Definition: VoigtArray.hpp:40
VoigtContainer(std::initializer_list< T > value)
Definition: VoigtArray.hpp:74
std::array< T, 6 > data_
Definition: VoigtArray.hpp:90
VoigtContainer()=default
static constexpr auto unique_indices
Definition: VoigtArray.hpp:54
VoigtContainer(const Array &array)
constexpr std::size_t size() const
Definition: VoigtArray.hpp:87
static constexpr auto indices
Definition: VoigtArray.hpp:42
const T & operator[](const VoigtIndex idx) const
Definition: VoigtArray.hpp:81
static constexpr auto diag_indices
Definition: VoigtArray.hpp:63
Definition: blackoilboundaryratevector.hh:39
VoigtIndex
Definition: VoigtArray.hpp:32