mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
added: voigt array
this is a container for vector data stored using voigt notation
This commit is contained in:
parent
a799f463fb
commit
5a3373f93f
@ -160,6 +160,7 @@ list (APPEND MAIN_SOURCE_FILES
|
||||
opm/simulators/utils/PressureAverage.cpp
|
||||
opm/simulators/utils/SerializationPackers.cpp
|
||||
opm/simulators/utils/UnsupportedFlowKeywords.cpp
|
||||
opm/simulators/utils/VoigtArray.cpp
|
||||
opm/simulators/utils/compressPartition.cpp
|
||||
opm/simulators/utils/gatherDeferredLogger.cpp
|
||||
opm/simulators/utils/phaseUsageFromDeck.cpp
|
||||
@ -419,6 +420,7 @@ list (APPEND TEST_SOURCE_FILES
|
||||
tests/test_stoppedwells.cpp
|
||||
tests/test_timer.cpp
|
||||
tests/test_vfpproperties.cpp
|
||||
tests/test_VoigtArray.cpp
|
||||
tests/test_wellmodel.cpp
|
||||
tests/test_wellprodindexcalculator.cpp
|
||||
tests/test_wellstate.cpp
|
||||
@ -985,6 +987,7 @@ list (APPEND PUBLIC_HEADER_FILES
|
||||
opm/simulators/utils/ParallelSerialization.hpp
|
||||
opm/simulators/utils/readDeck.hpp
|
||||
opm/simulators/utils/satfunc/RelpermDiagnostics.hpp
|
||||
opm/simulators/utils/VoigtArray.hpp
|
||||
opm/simulators/wells/ALQState.hpp
|
||||
opm/simulators/wells/BlackoilWellModel.hpp
|
||||
opm/simulators/wells/BlackoilWellModel_impl.hpp
|
||||
|
85
opm/simulators/utils/VoigtArray.cpp
Normal file
85
opm/simulators/utils/VoigtArray.cpp
Normal file
@ -0,0 +1,85 @@
|
||||
/*
|
||||
Copyright 2025 Equinor ASA
|
||||
|
||||
This file is part of the Open Porous Media Project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
#include <config.h>
|
||||
#include <opm/simulators/utils/VoigtArray.hpp>
|
||||
|
||||
#include <dune/common/fvector.hh>
|
||||
|
||||
#include <algorithm>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
template<class T>
|
||||
template<class Array>
|
||||
VoigtContainer<T>::
|
||||
VoigtContainer(const Array& array)
|
||||
{
|
||||
std::copy(array.begin(), array.end(), data_.begin());
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
VoigtArray<Scalar>::
|
||||
VoigtArray(const std::size_t size)
|
||||
{
|
||||
this->resize(size);
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void VoigtArray<Scalar>::
|
||||
resize(const std::size_t size)
|
||||
{
|
||||
std::for_each(this->data_.begin(), this->data_.end(),
|
||||
[size](auto& d) { d.resize(size); });
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void VoigtArray<Scalar>::
|
||||
assign(const std::size_t i, const VoigtContainer<Scalar>& array)
|
||||
{
|
||||
for (const auto idx : this->unique_indices) {
|
||||
(*this)[idx][i] = array[idx];
|
||||
}
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
Scalar VoigtArray<Scalar>::
|
||||
operator()(const VoigtIndex idx, const std::size_t i) const
|
||||
{
|
||||
return (*this)[idx].at(i);
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
Scalar& VoigtArray<Scalar>::
|
||||
operator()(const VoigtIndex idx, const std::size_t i)
|
||||
{
|
||||
return (*this)[idx].at(i);
|
||||
}
|
||||
|
||||
#define INSTANTIATE_TYPE(T) \
|
||||
template class VoigtArray<T>; \
|
||||
template VoigtContainer<T>::VoigtContainer(const std::array<T,6>&); \
|
||||
template VoigtContainer<T>::VoigtContainer(const Dune::FieldVector<T,6>&);
|
||||
|
||||
INSTANTIATE_TYPE(double)
|
||||
|
||||
#if FLOW_INSTANTIATE_FLOAT
|
||||
INSTANTIATE_TYPE(float)
|
||||
#endif
|
||||
|
||||
} // namespace Opm
|
112
opm/simulators/utils/VoigtArray.hpp
Normal file
112
opm/simulators/utils/VoigtArray.hpp
Normal file
@ -0,0 +1,112 @@
|
||||
/*
|
||||
Copyright 2025 Equinor ASA
|
||||
|
||||
This file is part of the Open Porous Media Project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#ifndef OPM_UTIL_VOIGT_ARRAY_HPP
|
||||
#define OPM_UTIL_VOIGT_ARRAY_HPP
|
||||
|
||||
#include <algorithm>
|
||||
#include <array>
|
||||
#include <cstddef>
|
||||
#include <initializer_list>
|
||||
#include <type_traits>
|
||||
#include <vector>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
enum class VoigtIndex {
|
||||
XX = 0, XY = 5, XZ = 4,
|
||||
YX = XY, YY = 1, YZ = 3,
|
||||
ZX = XZ, ZY = YZ, ZZ = 2,
|
||||
};
|
||||
|
||||
template<class T>
|
||||
class VoigtContainer
|
||||
{
|
||||
public:
|
||||
static constexpr auto indices = std::array{
|
||||
Opm::VoigtIndex::XX,
|
||||
Opm::VoigtIndex::XY,
|
||||
Opm::VoigtIndex::XZ,
|
||||
Opm::VoigtIndex::YX,
|
||||
Opm::VoigtIndex::YY,
|
||||
Opm::VoigtIndex::YZ,
|
||||
Opm::VoigtIndex::ZX,
|
||||
Opm::VoigtIndex::ZY,
|
||||
Opm::VoigtIndex::ZZ,
|
||||
};
|
||||
|
||||
static constexpr auto unique_indices = std::array{
|
||||
Opm::VoigtIndex::XX,
|
||||
Opm::VoigtIndex::YY,
|
||||
Opm::VoigtIndex::ZZ,
|
||||
Opm::VoigtIndex::YZ,
|
||||
Opm::VoigtIndex::XZ,
|
||||
Opm::VoigtIndex::XY
|
||||
};
|
||||
|
||||
static constexpr auto diag_indices = std::array{
|
||||
Opm::VoigtIndex::XX,
|
||||
Opm::VoigtIndex::YY,
|
||||
Opm::VoigtIndex::ZZ,
|
||||
};
|
||||
|
||||
VoigtContainer() = default;
|
||||
|
||||
template<class Array>
|
||||
VoigtContainer(const Array& array);
|
||||
|
||||
VoigtContainer(std::initializer_list<T> value)
|
||||
{
|
||||
std::copy_n(value.begin(),
|
||||
std::min(data_.size(), value.size()),
|
||||
data_.begin());
|
||||
}
|
||||
|
||||
const T& operator[](const VoigtIndex idx) const
|
||||
{ return data_[static_cast<std::underlying_type_t<VoigtIndex>>(idx)]; }
|
||||
|
||||
T& operator [](const VoigtIndex idx)
|
||||
{ return data_[static_cast<std::underlying_type_t<VoigtIndex>>(idx)]; }
|
||||
|
||||
constexpr std::size_t size() const { return data_.size(); }
|
||||
|
||||
protected:
|
||||
std::array<T, 6> data_{};
|
||||
};
|
||||
|
||||
template<class T> using SymmTensor = VoigtContainer<T>;
|
||||
|
||||
template<class Scalar>
|
||||
class VoigtArray : public VoigtContainer<std::vector<Scalar>>
|
||||
{
|
||||
public:
|
||||
VoigtArray() = default;
|
||||
explicit VoigtArray(const std::size_t size);
|
||||
|
||||
void resize(const std::size_t size);
|
||||
|
||||
Scalar operator()(const VoigtIndex idx, const std::size_t i) const;
|
||||
Scalar& operator()(const VoigtIndex idx, const std::size_t i);
|
||||
|
||||
void assign(const std::size_t i, const VoigtContainer<Scalar>& array);
|
||||
};
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
#endif // OPM_UTIL_VOIGT_ARRAY_HPP
|
87
tests/test_VoigtArray.cpp
Normal file
87
tests/test_VoigtArray.cpp
Normal file
@ -0,0 +1,87 @@
|
||||
/*
|
||||
Copyright 2025 Equinor ASA.
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <config.h>
|
||||
|
||||
#include <opm/simulators/utils/VoigtArray.hpp>
|
||||
|
||||
#define BOOST_TEST_MODULE VoigArrayTest
|
||||
#include <boost/mpl/list.hpp>
|
||||
#include <boost/test/unit_test.hpp>
|
||||
|
||||
#include <stdexcept>
|
||||
|
||||
namespace {
|
||||
#if FLOW_INSTANTIATE_FLOAT
|
||||
using Types = boost::mpl::list<float,double>;
|
||||
#else
|
||||
using Types = boost::mpl::list<double>;
|
||||
#endif
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE_TEMPLATE(DefaultConstructed, Scalar, Types)
|
||||
{
|
||||
Opm::VoigtArray<Scalar> array;
|
||||
static constexpr auto& indices = Opm::VoigtArray<Scalar>::indices;
|
||||
BOOST_CHECK_EQUAL(array.size(), 6);
|
||||
for (const auto i : indices) {
|
||||
BOOST_CHECK(array[i].empty());
|
||||
}
|
||||
|
||||
for (const auto i : indices) {
|
||||
BOOST_CHECK(array[i].empty());
|
||||
BOOST_CHECK_THROW(array(i, 1), std::out_of_range);
|
||||
}
|
||||
|
||||
array.resize(5);
|
||||
for (const auto i : indices) {
|
||||
BOOST_CHECK(!array[i].empty());
|
||||
BOOST_CHECK_NO_THROW(array(i, 1));
|
||||
}
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE_TEMPLATE(SizeConstructed, Scalar, Types)
|
||||
{
|
||||
Opm::VoigtArray<Scalar> array(4);
|
||||
BOOST_CHECK_EQUAL(array.size(), 6);
|
||||
static constexpr auto& indices = Opm::VoigtArray<Scalar>::indices;
|
||||
static constexpr auto& uindices = Opm::VoigtArray<Scalar>::unique_indices;
|
||||
for (const auto i : indices) {
|
||||
BOOST_CHECK_EQUAL(array[i].size(), 4);
|
||||
}
|
||||
|
||||
for (const auto i : uindices) {
|
||||
std::fill(array[i].begin(), array[i].end(), static_cast<Scalar>(i));
|
||||
}
|
||||
|
||||
for (const auto i: indices) {
|
||||
BOOST_CHECK_EQUAL(array[i].size(), 4);
|
||||
BOOST_CHECK_EQUAL(array[i][1], static_cast<Scalar>(i));
|
||||
BOOST_CHECK_EQUAL(array(i,1), static_cast<Scalar>(i));
|
||||
}
|
||||
|
||||
BOOST_CHECK_EQUAL(array(Opm::VoigtIndex::XY, 1),
|
||||
array(Opm::VoigtIndex::YX, 1));
|
||||
|
||||
BOOST_CHECK_EQUAL(array(Opm::VoigtIndex::XZ, 2),
|
||||
array(Opm::VoigtIndex::ZX, 2));
|
||||
|
||||
BOOST_CHECK_EQUAL(array(Opm::VoigtIndex::YZ, 3),
|
||||
array(Opm::VoigtIndex::ZY, 3));
|
||||
}
|
Loading…
Reference in New Issue
Block a user