preCICE
Loading...
Searching...
No Matches
OnBoundViolationTest.cpp
Go to the documentation of this file.
1#ifndef PRECICE_NO_MPI
2
3#include <Eigen/Core>
4#include <Eigen/Dense>
5
10#include "testing/Meshes.hpp"
12#include "testing/Testing.hpp"
13
14BOOST_AUTO_TEST_SUITE(AccelerationTests)
15
16using namespace precice;
17using namespace precice::acceleration;
19
21BOOST_AUTO_TEST_CASE(CheckBoundViolation)
22{
24
25 const int scalarDataDimension = 1;
26 const int vectorDataDimension = 3;
27
28 auto dummyMesh = testing::makeDummy2DMesh(2);
29
30 // only set lower bound for scalar data
31 std::vector<std::optional<double>> lowerBoundedL = {0};
32 std::vector<std::optional<double>> lowerBoundedU(scalarDataDimension);
33 // only set upper bound for scalar data
34 std::vector<std::optional<double>> upperBoundedL(scalarDataDimension);
35 std::vector<std::optional<double>> upperBoundedU = {0};
36 // set both lower and upper bound for scalar data
37 std::vector<std::optional<double>> lowerUpperBoundedL = {0};
38 std::vector<std::optional<double>> lowerUpperBoundedU = {1};
39 // set upper bound for z component in vector data
40 std::vector<std::optional<double>> zUpperBoundedL(vectorDataDimension);
41 std::vector<std::optional<double>> zUpperBoundedU(vectorDataDimension);
42 zUpperBoundedU[2] = 1;
43 // set lower bound for x-component
44 std::vector<std::optional<double>> xLowerBoundedL(vectorDataDimension);
45 std::vector<std::optional<double>> xLowerBoundedU(vectorDataDimension);
46 xLowerBoundedL[0] = 0;
47 // set upper bound for y-component
48 std::vector<std::optional<double>> yUpperBoundedL(vectorDataDimension);
49 std::vector<std::optional<double>> yUpperBoundedU(vectorDataDimension);
50 yUpperBoundedU[1] = 1;
51 // set upper bound for x-component
52 std::vector<std::optional<double>> xUpperBoundedL(vectorDataDimension);
53 std::vector<std::optional<double>> xUpperBoundedU(vectorDataDimension);
54 xUpperBoundedU[0] = 1;
55
56 mesh::PtrData lowerBounded(new mesh::Data("lowerBounded", 0, scalarDataDimension, 1, 2, lowerBoundedL, lowerBoundedU));
57 lowerBounded->emplaceSampleAtTime(0.0, {0., 0.});
58 mesh::PtrData upperBounded(new mesh::Data("upperBounded", 1, scalarDataDimension, 1, 2, upperBoundedL, upperBoundedU));
59 upperBounded->emplaceSampleAtTime(0.0, {0., 0.});
60 mesh::PtrData lowerUpperBounded(new mesh::Data("lowerUpperBounded", 2, scalarDataDimension, 1, 2, lowerUpperBoundedL, lowerUpperBoundedU));
61 lowerUpperBounded->emplaceSampleAtTime(0.0, {0., 0.});
62
63 mesh::PtrData zUpperBounded(new mesh::Data("zUpperBounded", 3, vectorDataDimension, 1, 2, zUpperBoundedL, zUpperBoundedU));
64 zUpperBounded->emplaceSampleAtTime(0.0, {0., 0., 0., 0., 0., 0.});
65 mesh::PtrData xLowerBounded(new mesh::Data("xLowerBounded", 4, vectorDataDimension, 1, 2, xLowerBoundedL, xLowerBoundedU));
66 xLowerBounded->emplaceSampleAtTime(0.0, {0., 0., 0., 0., 0., 0.});
67 mesh::PtrData yUpperBounded(new mesh::Data("yUpperBounded", 5, vectorDataDimension, 1, 2, yUpperBoundedL, yUpperBoundedU));
68 yUpperBounded->emplaceSampleAtTime(0.0, {0., 0., 0., 0., 0., 0.});
69 mesh::PtrData xUpperBounded(new mesh::Data("xUpperBounded", 6, vectorDataDimension, 1, 2, xUpperBoundedL, xUpperBoundedU));
70 xUpperBounded->emplaceSampleAtTime(0.0, {0., 0., 0., 0., 0., 0.});
71
72 cplscheme::PtrCouplingData lowerBoundedPtr = makeCouplingData(lowerBounded, dummyMesh, true);
73 cplscheme::PtrCouplingData upperBoundedPtr = makeCouplingData(upperBounded, dummyMesh, true);
74 cplscheme::PtrCouplingData lowerUpperBoundedPtr = makeCouplingData(lowerUpperBounded, dummyMesh, true);
75 cplscheme::PtrCouplingData zUpperBoundedPtr = makeCouplingData(zUpperBounded, dummyMesh, true);
76 cplscheme::PtrCouplingData xLowerBoundedPtr = makeCouplingData(xLowerBounded, dummyMesh, true);
77 cplscheme::PtrCouplingData yUpperBoundedPtr = makeCouplingData(yUpperBounded, dummyMesh, true);
78 cplscheme::PtrCouplingData xUpperBoundedPtr = makeCouplingData(xUpperBounded, dummyMesh, true);
79 for (auto &cd : {lowerBoundedPtr, upperBoundedPtr, lowerUpperBoundedPtr, zUpperBoundedPtr, xLowerBoundedPtr, yUpperBoundedPtr, xUpperBoundedPtr}) {
80 cd->storeIteration();
81 }
82
83 cplscheme::DataMap cplData{
84 {0, lowerBoundedPtr},
85 {1, upperBoundedPtr},
86 {2, lowerUpperBoundedPtr},
87 {3, zUpperBoundedPtr},
88 {4, xLowerBoundedPtr},
89 {5, yUpperBoundedPtr},
90 {6, xUpperBoundedPtr}};
91
92 IQNILSAcceleration acceleration(1.0, false, 10, 0, 0, 1e-10, {0}, acceleration::Acceleration::OnBoundViolation::Clamp, nullptr, false);
93 acceleration.initialize(cplData);
94 Eigen::VectorXd testData(3 * 2 + 4 * 6);
95 testData << 2.0, -0.1, 2.0, 0.1, -0.1, 1.1, 0.1, 0.1, 1.1, 0.1, 0.1, 0.1, 0.1, 0.1, 0.1, -0.1, 0.1, 0.1, 0.1, 0.1, 0.1, 0.1, 1.1, 0.1, 1.1, 0.1, 0.1, 0.1, 0.1, 0.1;
96 auto violatingIDs = acceleration.checkBoundViolation(testData, cplData);
97
98 std::vector<int> expectedViolations{0, 1, 2, 3, 4, 5, 6};
99 BOOST_TEST(violatingIDs == expectedViolations, boost::test_tools::per_element());
100}
101
102PRECICE_TEST_SETUP(1_rank)
103BOOST_AUTO_TEST_CASE(testOnBoundViolationClamp)
104{
105 PRECICE_TEST();
106
107 Eigen::VectorXd _data(6);
108 _data << 0.1, 0.3, 0.4, -0.4, 1.2, 0.5;
109 Eigen::Map<Eigen::MatrixXd> dataMap(_data.data(), 2, 3);
110
111 Eigen::VectorXd _dataAfterUpdate(6);
112 _dataAfterUpdate << 0.1, 0.3, 0.4, 0.0, 0.5, 0.5;
113 Eigen::Map<Eigen::MatrixXd> dataAfterUpdateMap(_dataAfterUpdate.data(), 2, 3);
114
115 std::vector<std::optional<double>> _lowerBound(2);
116 _lowerBound[0] = 0;
117 _lowerBound[1] = 0;
118 std::vector<std::optional<double>> _upperBound(2);
119 _upperBound[0] = 0.5;
120 _upperBound[1] = 1;
121
122 IQNILSAcceleration acceleration(1.0, false, 10, 0, 0, 1e-10, {0}, acceleration::Acceleration::OnBoundViolation::Clamp, nullptr, false);
123 acceleration.clampToBounds(dataMap, _lowerBound, _upperBound);
124
125 // check if the violating values have been correctly truncated
126 BOOST_TEST(testing::equals(dataMap, dataAfterUpdateMap, 1e-10));
127}
128
129PRECICE_TEST_SETUP(1_rank)
130BOOST_AUTO_TEST_CASE(testOnBoundViolationScale)
131{
132 PRECICE_TEST();
133
134 Eigen::VectorXd _data(6);
135 _data << 0.1, 0.3, 0.4, -0.4, 1.2, 0.5;
136 Eigen::Map<Eigen::MatrixXd> dataMap(_data.data(), 2, 3);
137
138 Eigen::VectorXd _dataUpdate(6);
139 _dataUpdate << -0.6, -0.3, 0.3, -0.6, 0.8, 0.12;
140 Eigen::Map<Eigen::MatrixXd> dataUpdateMap(_dataUpdate.data(), 2, 3);
141
142 std::vector<std::optional<double>> _lowerBound(2);
143 _lowerBound[0] = 0;
144 _lowerBound[1] = 0;
145 std::vector<std::optional<double>> _upperBound(2);
146 _upperBound[0] = 0.5;
147 _upperBound[1] = 1;
148
149 IQNILSAcceleration acceleration(1.0, false, 10, 0, 0, 1e-10, {0}, acceleration::Acceleration::OnBoundViolation::ScaleToBound, nullptr, false);
150 auto scaleFactor = acceleration.computeShorteningFactor(dataMap, dataUpdateMap, _lowerBound, _upperBound);
151 BOOST_TEST(testing::equals(scaleFactor, 0.875, 1e-10));
152}
153
154PRECICE_TEST_SETUP(1_rank)
155BOOST_AUTO_TEST_CASE(testOnBoundViolationScaleDivZero)
156{
157 PRECICE_TEST();
158
159 Eigen::VectorXd _data(6);
160 _data << -1e-13, -0.3, 1e-13, -1e-13, 1e-13, 1e-13;
161 Eigen::Map<Eigen::MatrixXd> dataMap(_data.data(), 2, 3);
162
163 Eigen::VectorXd _dataUpdate(6);
164 _dataUpdate << -2e-13, -0.3, 0, 3e-13, -2e-13, 2e-13;
165 Eigen::Map<Eigen::MatrixXd> dataUpdateMap(_dataUpdate.data(), 2, 3);
166
167 std::vector<std::optional<double>> _lowerBound(2);
168 _lowerBound[0] = 0;
169 _lowerBound[1] = -1;
170 std::vector<std::optional<double>> _upperBound(2);
171 _upperBound[0] = 0.5;
172 _upperBound[1] = 0;
173
174 IQNILSAcceleration acceleration(1.0, false, 10, 0, 0, 1e-10, {0}, acceleration::Acceleration::OnBoundViolation::ScaleToBound, nullptr, false);
175 auto scaleFactor = acceleration.computeShorteningFactor(dataMap, dataUpdateMap, _lowerBound, _upperBound);
176 BOOST_TEST(testing::equals(scaleFactor, 0.5, 1e-10));
177}
178
180
181#endif
BOOST_AUTO_TEST_CASE(testIQNIMVJPPWithSubsteps)
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
cplscheme::PtrCouplingData makeCouplingData(mesh::PtrData data, mesh::PtrMesh mesh, bool exchangeSubsteps)
Definition helper.hpp:10
#define PRECICE_TEST()
Definition Testing.hpp:39
#define PRECICE_TEST_SETUP(...)
Creates and attaches a TestSetup to a Boost test case.
Definition Testing.hpp:29
Interface quasi-Newton with interface least-squares approximation.
Describes a set of data values belonging to the vertices of a mesh.
Definition Data.hpp:26
contains implementations of acceleration schemes.
std::shared_ptr< CouplingData > PtrCouplingData
std::map< int, PtrCouplingData > DataMap
std::shared_ptr< Data > PtrData
boost::test_tools::predicate_result equals(const std::vector< float > &VectorA, const std::vector< float > &VectorB, float tolerance)
equals to be used in tests. Compares two std::vectors using a given tolerance. Prints both operands o...
Definition Testing.cpp:93
cplscheme::PtrCouplingData makeCouplingData(mesh::PtrData data, mesh::PtrMesh mesh, bool exchangeSubsteps)
Definition helper.hpp:10
auto makeDummy2DMesh(size_t nVertices)
Definition Meshes.hpp:21
Main namespace of the precice library.