-
Notifications
You must be signed in to change notification settings - Fork 9
/
test_integration_da.cc
243 lines (205 loc) · 8.14 KB
/
test_integration_da.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
/* Copyright (c) 2016 - 2024, the adamantine authors.
*
* This file is subject to the Modified BSD License and may not be distributed
* without copyright and license information. Please refer to the file LICENSE
* for the text and further information on this license.
*/
#define BOOST_TEST_MODULE Integration_Data_Assimilation
#include "../application/adamantine.hh"
#include <boost/property_tree/info_parser.hpp>
#include <filesystem>
#include <fstream>
#include "main.cc"
namespace tt = boost::test_tools;
double integration_da(MPI_Comm communicator, bool l2_norm)
{
std::vector<adamantine::Timer> timers;
initialize_timers(communicator, timers);
// Read the input.
std::string const filename = "bare_plate_L_da.info";
adamantine::ASSERT_THROW(std::filesystem::exists(filename) == true,
"The file " + filename + " does not exist.");
boost::property_tree::ptree database;
boost::property_tree::info_parser::read_info(filename, database);
// Run the simulation
auto result = run_ensemble<3, dealii::MemorySpace::Host>(communicator,
database, timers);
if (l2_norm)
{
double norm = -1.;
for (auto const &r : result)
{
norm = std::max(r.l2_norm(), norm);
}
return dealii::Utilities::MPI::max(norm, communicator);
}
// Three ensemble members expected
unsigned int local_result_size = result.size();
unsigned int global_result_size =
dealii::Utilities::MPI::sum(local_result_size, communicator);
BOOST_TEST(global_result_size == 3);
// Get the average minimum value for each ensemble member, which is very close
// to the initial temperature
double sum = 0.0;
for (unsigned int member = 0; member < local_result_size; ++member)
{
double min_val = std::numeric_limits<double>::max();
for (unsigned int i = 0;
i < result.at(member).block(0).locally_owned_size(); ++i)
{
if (result.at(member).block(0).local_element(i) < min_val)
min_val = result.at(member).block(0).local_element(i);
}
sum += min_val;
}
double partial_average_minimum_value = sum / global_result_size;
double average_minimum_value =
dealii::Utilities::MPI::sum(partial_average_minimum_value, communicator);
// Based on the experimental data, the expected temperature is ~200.0
BOOST_TEST(average_minimum_value >= 200.0);
BOOST_TEST(average_minimum_value < 300.0);
MPI_Barrier(communicator);
return average_minimum_value;
}
BOOST_AUTO_TEST_CASE(integration_data_assimilation)
{
bool l2_norm = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) > 3
? true
: false;
double result_world = integration_da(MPI_COMM_WORLD, l2_norm);
double result_self = integration_da(MPI_COMM_SELF, l2_norm);
BOOST_TEST(result_world == result_self, tt::tolerance(1e-12));
}
double integration_da_point_cloud_add_mat(MPI_Comm communicator, bool l2_norm)
{
std::vector<adamantine::Timer> timers;
initialize_timers(communicator, timers);
// Read the input.
std::string const filename = "integration_da_add_material.info";
adamantine::ASSERT_THROW(std::filesystem::exists(filename) == true,
"The file " + filename + " does not exist.");
boost::property_tree::ptree database;
boost::property_tree::info_parser::read_info(filename, database);
// Run the simulation
auto result = run_ensemble<3, dealii::MemorySpace::Host>(communicator,
database, timers);
if (l2_norm)
{
double norm = -1.;
for (auto const &r : result)
{
norm = std::max(r.l2_norm(), norm);
}
return dealii::Utilities::MPI::max(norm, communicator);
}
// Three ensemble members expected
unsigned int local_result_size = result.size();
unsigned int global_result_size =
dealii::Utilities::MPI::sum(local_result_size, communicator);
BOOST_TEST(global_result_size == 3);
// Get the average minimum value for each ensemble member, which is very close
// to the initial temperature
double sum = 0.0;
for (unsigned int member = 0; member < local_result_size; ++member)
{
double min_val = std::numeric_limits<double>::max();
for (unsigned int i = 0;
i < result.at(member).block(0).locally_owned_size(); ++i)
{
if (result.at(member).block(0).local_element(i) < min_val)
min_val = result.at(member).block(0).local_element(i);
}
sum += min_val;
}
double partial_average_minimum_value = sum / global_result_size;
double average_minimum_value =
dealii::Utilities::MPI::sum(partial_average_minimum_value, communicator);
// Based on the experimental data, the expected temperature is ~200.0
BOOST_CHECK(average_minimum_value >= 200.0);
BOOST_CHECK(average_minimum_value < 300.0);
MPI_Barrier(communicator);
return average_minimum_value;
}
BOOST_AUTO_TEST_CASE(integration_3D_da_point_cloud_add_material)
{
/*
* This integration test checks that the data assimilation using point cloud
* data works while adding material.
*/
bool l2_norm = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) > 3
? true
: false;
double result_world =
integration_da_point_cloud_add_mat(MPI_COMM_WORLD, l2_norm);
double result_self =
integration_da_point_cloud_add_mat(MPI_COMM_SELF, l2_norm);
BOOST_TEST(result_world == result_self, tt::tolerance(1e-12));
}
double integration_da_ray_add_mat(MPI_Comm communicator, bool l2_norm)
{
std::vector<adamantine::Timer> timers;
initialize_timers(communicator, timers);
// Read the input.
std::string const filename = "integration_da_add_material.info";
adamantine::ASSERT_THROW(std::filesystem::exists(filename) == true,
"The file " + filename + " does not exist.");
boost::property_tree::ptree database;
boost::property_tree::info_parser::read_info(filename, database);
// Change experiment file name
database.put("experiment.file",
"integration_da_add_material_expt_ray_#camera_#frame.csv");
// Change experiment data format to ray
database.put("experiment.format", "ray");
// Run the simulation
auto result = run_ensemble<3, dealii::MemorySpace::Host>(communicator,
database, timers);
if (l2_norm)
{
double norm = -1.;
for (auto const &r : result)
{
norm = std::max(r.l2_norm(), norm);
}
return dealii::Utilities::MPI::max(norm, communicator);
}
// Three ensemble members expected
unsigned int local_result_size = result.size();
unsigned int global_result_size =
dealii::Utilities::MPI::sum(local_result_size, communicator);
BOOST_TEST(global_result_size == 3);
// Get the average minimum value for each ensemble member, which is very close
// to the initial temperature
double sum = 0.0;
for (unsigned int member = 0; member < result.size(); ++member)
{
double min_val = std::numeric_limits<double>::max();
for (unsigned int i = 0;
i < result.at(member).block(0).locally_owned_size(); ++i)
{
if (result.at(member).block(0).local_element(i) < min_val)
min_val = result.at(member).block(0).local_element(i);
}
sum += min_val;
}
double partial_average_minimum_value = sum / global_result_size;
double average_minimum_value =
dealii::Utilities::MPI::sum(partial_average_minimum_value, communicator);
// Based on the experimental data, the expected temperature is ~200.0
BOOST_CHECK(average_minimum_value >= 200.0);
BOOST_CHECK(average_minimum_value < 300.0);
MPI_Barrier(communicator);
return average_minimum_value;
}
BOOST_AUTO_TEST_CASE(integration_3D_da_ray_add_material)
{
/*
* This integration test checks that the data assimilation using point cloud
* data works while adding material.
*/
bool l2_norm = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) > 3
? true
: false;
double result_world = integration_da_ray_add_mat(MPI_COMM_WORLD, l2_norm);
double result_self = integration_da_ray_add_mat(MPI_COMM_SELF, l2_norm);
BOOST_TEST(result_world == result_self, tt::tolerance(1e-12));
}