-
Notifications
You must be signed in to change notification settings - Fork 7
/
test_integration.cpp
540 lines (475 loc) · 21.5 KB
/
test_integration.cpp
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
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
#define BOOST_TEST_DYN_LINK
// to resolve https://github.com/open-mpi/ompi/issues/5157
#define OMPI_SKIP_MPICXX 1
#include <mpi.h>
#include <boost/serialization/export.hpp>
#include <boost/test/unit_test.hpp>
#include <cstdio>
#include "TaskCount.hpp"
#include "combischeme/CombiMinMaxScheme.hpp"
#include "io/H5InputOutput.hpp"
#include "loadmodel/LearningLoadModel.hpp"
#include "loadmodel/LinearLoadModel.hpp"
#include "manager/CombiParameters.hpp"
#include "manager/ProcessGroupManager.hpp"
#include "manager/ProcessGroupWorker.hpp"
#include "manager/ProcessManager.hpp"
#include "sparsegrid/DistributedSparseGridUniform.hpp"
#include "task/Task.hpp"
#include "utils/Config.hpp"
#include "utils/MonteCarlo.hpp"
#include "utils/Types.hpp"
#include "stdlib.h"
#include "test_helper.hpp"
using namespace combigrid;
// omitted here, because already defined in test_thirdLevel.cpp
// BOOST_CLASS_EXPORT(TaskCount)
bool checkReducedFullGridIntegration(ProcessGroupWorker& worker, int nrun) {
const auto& tasks = worker.getTasks();
int numGrids = (int)worker.getCombiParameters().getNumGrids();
BOOST_CHECK(tasks.size() > 0);
BOOST_CHECK(numGrids > 0);
// to check if any data was actually compared
bool any = false;
for (const auto& t : tasks) {
for (int g = 0; g < numGrids; g++) {
DistributedFullGrid<CombiDataType>& dfg = t->getDistributedFullGrid(g);
for (auto b : dfg.returnBoundaryFlags()) {
BOOST_CHECK(b == 2);
}
TestFnCount<CombiDataType> initialFunction;
for (IndexType li = 0; li < dfg.getNrLocalElements(); ++li) {
std::vector<double> coords(dfg.getDimension());
dfg.getCoordsLocal(li, coords);
CombiDataType expected = initialFunction(coords, nrun);
CombiDataType occuring = dfg.getData()[li];
for (auto& c : coords) {
BOOST_CHECK(c >= 0.);
BOOST_CHECK(c <= 1.);
}
// checking absolute and real value since comparing std::complex may be tricky
BOOST_CHECK_CLOSE(std::abs(expected), std::abs(occuring), TestHelper::tolerance);
BOOST_CHECK_CLOSE(std::real(expected), std::real(occuring), TestHelper::tolerance);
// BOOST_REQUIRE_CLOSE(expected, occuring, TestHelper::tolerance);
any = true;
}
}
}
BOOST_REQUIRE(any);
return any;
}
void checkIntegration(size_t ngroup = 1, size_t nprocs = 1, BoundaryType boundaryV = 2,
bool pretendThirdLevel = true) {
size_t size = ngroup * nprocs + 1;
BOOST_REQUIRE(TestHelper::checkNumMPIProcsAvailable(size));
CommunicatorType comm = TestHelper::getComm(size);
if (comm == MPI_COMM_NULL) {
BOOST_TEST_CHECKPOINT("drop out of test comm");
return;
}
BOOST_TEST_CHECKPOINT("initialize stats");
combigrid::Stats::initialize();
theMPISystem()->initWorldReusable(comm, ngroup, nprocs);
// theMPISystem()->init(ngroup, nprocs);
DimType dim = 2;
LevelVector lmin(dim, 2);
LevelVector lmax(dim, 5);
size_t ncombi = 4;
BOOST_CHECK_EQUAL(theMPISystem()->getWorldSize(), size);
BOOST_CHECK_EQUAL(getCommSize(theMPISystem()->getWorldComm()), size);
WORLD_MANAGER_EXCLUSIVE_SECTION {
// make sure the manager's ranks are set right
BOOST_CHECK_EQUAL(getCommRank(theMPISystem()->getWorldComm()), size - 1);
BOOST_CHECK_EQUAL(getCommSize(theMPISystem()->getGlobalComm()), ngroup + 1);
BOOST_CHECK_EQUAL(getCommRank(theMPISystem()->getGlobalComm()), ngroup);
ProcessGroupManagerContainer pgroups;
for (int i = 0; i < ngroup; ++i) {
int pgroupRootID(i);
pgroups.emplace_back(std::make_shared<ProcessGroupManager>(pgroupRootID));
}
auto loadmodel = std::unique_ptr<LoadModel>(new LinearLoadModel());
std::vector<BoundaryType> boundary(dim, boundaryV);
CombiMinMaxScheme combischeme(dim, lmin, lmax);
combischeme.createAdaptiveCombischeme();
std::vector<LevelVector> levels = combischeme.getCombiSpaces();
std::vector<combigrid::real> coeffs = combischeme.getCoeffs();
// create Tasks
TaskContainer tasks;
std::vector<size_t> taskIDs;
for (size_t i = 0; i < levels.size(); i++) {
Task* t = new TaskCount(levels[i], boundary, coeffs[i], loadmodel.get());
tasks.push_back(t);
taskIDs.push_back(t->getID());
}
// why are there duplicate task IDs over different calls to this function??
// ah its because different ranks will be master, so there are different Task::count instances
// std::cout << "taskIDs " << taskIDs << std::endl;
// create combiparameters
BOOST_TEST_CHECKPOINT("manager create combi parameters");
CombiParameters params(dim, lmin, lmax, boundary, levels, coeffs, taskIDs, ncombi, 1,
CombinationVariant::sparseGridReduce, {static_cast<int>(nprocs), 1},
LevelVector(0), LevelVector(0), 16, false);
if (nprocs == 5 && boundaryV == 2) {
params.setDecomposition({{0, 6, 13, 20, 27}, {0}});
} else if (nprocs == 4 && boundaryV == 2) {
// should be the same as default decomposition with forwardDecomposition
params.setDecomposition({{0, 9, 17, 25}, {0}});
} else if (nprocs == 3) {
params.setDecomposition({{0, 15, 20}, {0}});
}
// create abstraction for Manager
ProcessManager manager{pgroups, tasks, params, std::move(loadmodel)};
// the combiparameters are sent to all process groups before the
// computations start
BOOST_TEST_CHECKPOINT("manager update combi parameters");
manager.updateCombiParameters();
/* distribute task according to load model and start computation for
* the first time */
BOOST_TEST_CHECKPOINT("run first");
auto start = std::chrono::high_resolution_clock::now();
if (pretendThirdLevel) {
manager.runfirst(true);
manager.pretendUnifySubspaceSizesThirdLevel();
} else {
manager.runfirst();
}
auto end = std::chrono::high_resolution_clock::now();
auto duration = std::chrono::duration_cast<std::chrono::milliseconds>(end - start);
BOOST_TEST_MESSAGE("manager run first solver step: " << duration.count() << " milliseconds");
for (size_t it = 0; it < ncombi - 1; ++it) {
BOOST_TEST_CHECKPOINT("combine");
start = std::chrono::high_resolution_clock::now();
manager.combine();
if (pretendThirdLevel) {
manager.pretendCombineThirdLevelForWorkers();
}
end = std::chrono::high_resolution_clock::now();
duration = std::chrono::duration_cast<std::chrono::milliseconds>(end - start);
BOOST_TEST_MESSAGE("manager combine: " << duration.count() << " milliseconds");
BOOST_TEST_CHECKPOINT("run next");
start = std::chrono::high_resolution_clock::now();
manager.runnext();
end = std::chrono::high_resolution_clock::now();
duration = std::chrono::duration_cast<std::chrono::milliseconds>(end - start);
BOOST_TEST_MESSAGE("manager run: " << duration.count() << " milliseconds");
}
manager.combine();
Stats::startEvent("manager get norms");
// get all kinds of norms
manager.getLpNorms(0);
manager.getLpNorms(1);
manager.getLpNorms(2);
Stats::stopEvent("manager get norms");
BOOST_TEST_CHECKPOINT("write solution");
std::string filename("integration_" + std::to_string(ncombi) + ".raw");
BOOST_TEST_CHECKPOINT("write solution " + filename);
Stats::startEvent("manager write solution");
manager.parallelEval(lmax, filename, 0);
BOOST_TEST_CHECKPOINT("write min/max coefficients");
manager.writeSparseGridMinMaxCoefficients("integration_" + std::to_string(boundaryV) +
"_sparse_minmax");
Stats::stopEvent("manager write solution");
BOOST_TEST_MESSAGE("manager write solution: " << Stats::getDuration("manager write solution")
<< " milliseconds");
filename = "integration_" + std::to_string(boundaryV) + "_dsgs";
BOOST_TEST_CHECKPOINT("write DSGS " + filename);
Stats::startEvent("manager write DSG");
manager.writeDSGsToDisk(filename);
manager.readDSGsFromDisk(filename);
Stats::stopEvent("manager write DSG");
BOOST_TEST_MESSAGE("manager write/read DSG: " << Stats::getDuration("manager write DSG")
<< " milliseconds");
// test Monte-Carlo interpolation
// only if boundary values are used
if (boundaryV > 0) {
BOOST_TEST_CHECKPOINT("MC interpolation coordinates");
auto interpolationCoords = montecarlo::getRandomCoordinates(1000, dim);
BOOST_TEST_CHECKPOINT("MC interpolation");
Stats::startEvent("manager interpolate");
auto values = manager.interpolateValues(interpolationCoords);
Stats::stopEvent("manager interpolate");
BOOST_TEST_MESSAGE("manager interpolate: " << Stats::getDuration("manager interpolate")
<< " milliseconds");
if (boundaryV > 1) {
TestFnCount<CombiDataType> initialFunction;
for (size_t i = 0; i < interpolationCoords.size(); ++i) {
if (std::abs(initialFunction(interpolationCoords[i], ncombi) - values[i]) >
TestHelper::tolerance) {
std::cout << "err " << interpolationCoords.size() << interpolationCoords[i] << " " << i
<< std::endl;
}
auto ref = initialFunction(interpolationCoords[i], ncombi);
BOOST_CHECK_CLOSE(std::abs(ref), std::abs(values[i]), TestHelper::tolerance);
BOOST_CHECK_CLOSE(std::real(ref), std::real(values[i]), TestHelper::tolerance);
}
}
#ifdef DISCOTEC_USE_HIGHFIVE
// output files are not needed, remove them right away
// (if this doesn't happen, there may be hdf5 errors due to duplicate task IDs)
auto status = system("rm integration_interpolated*.h5");
BOOST_WARN_GE(status, 0);
sleep(1);
manager.writeInterpolatedValuesSingleFile(interpolationCoords, "integration_interpolated");
BOOST_TEST_CHECKPOINT("wrote interpolated values to single file");
manager.writeInterpolatedValuesPerGrid(interpolationCoords, "integration_interpolated");
BOOST_TEST_CHECKPOINT("wrote interpolated values per grid");
manager.writeInterpolationCoordinates(interpolationCoords, "integration_interpolated");
BOOST_TEST_CHECKPOINT("wrote interpolation coordinates");
sleep(1); // wait for filesystem to catch up
decltype(values) valuesAllGridsRead;
// the number of combinations that the worker counts is the number of calls to combine() and
// pretendCombineThirdLevelForWorkers()
h5io::readH5Values(valuesAllGridsRead, "integration_interpolated_values_" +
std::to_string(ncombi + (ncombi - 1)) + ".h5");
decltype(interpolationCoords) interpolationCoordsRead =
broadcastParameters::getCoordinatesFromRankZero("integration_interpolated_coords.h5",
MPI_COMM_SELF);
BOOST_TEST_CHECKPOINT("read interpolation coordinates");
BOOST_CHECK_EQUAL_COLLECTIONS(values.begin(), values.end(), valuesAllGridsRead.begin(),
valuesAllGridsRead.end());
for (size_t i = 0; i < interpolationCoords.size(); ++i) {
for (size_t d = 0; d < dim; ++d) {
BOOST_CHECK_EQUAL(interpolationCoords[i][d], interpolationCoordsRead[i][d]);
}
}
#endif // def DISCOTEC_USE_HIGHFIVE
}
manager.exit();
// if output files are not needed, remove them right away
remove(("integration_" + std::to_string(ncombi) + "_0.raw").c_str());
remove(("integration_" + std::to_string(ncombi) + "_0.raw_header").c_str());
BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getGlobalComm()));
}
else {
BOOST_CHECK_EQUAL(getCommSize(theMPISystem()->getLocalComm()), nprocs);
if (nprocs == 1) {
BOOST_CHECK(theMPISystem()->isProcessGroupMaster());
}
BOOST_TEST_CHECKPOINT("Worker starts");
BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getLocalComm()));
ProcessGroupWorker pgroup;
SignalType signal = -1;
// omitting to count RUN_FIRST signal, as it is executed once for every task
int nrun = 1;
while (signal != EXIT) {
BOOST_TEST_CHECKPOINT("Last Successful Worker Signal " + std::to_string(signal));
signal = pgroup.wait();
if (signal == RUN_NEXT) {
BOOST_CHECK(pgroup.getTasks().size() > 0);
++nrun;
}
if (signal == COMBINE) {
// after combination check workers' grids
// only if boundary values are used
if (boundaryV == 2) {
BOOST_CHECK(checkReducedFullGridIntegration(pgroup, nrun));
}
// write partial stats
if (theMPISystem()->getWorldRank() < nprocs) {
Stats::writePartial("integration_partial_timers_group.json",
theMPISystem()->getLocalComm());
}
}
}
BOOST_CHECK_EQUAL(nrun, ncombi);
BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getLocalComm()));
MASTER_EXCLUSIVE_SECTION { BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getGlobalComm())); }
}
combigrid::Stats::finalize();
Stats::write("integration_" + std::to_string(ngroup) + "_" + std::to_string(nprocs) + ".json");
MPI_Barrier(comm);
BOOST_CHECK(!TestHelper::testStrayMessages(comm));
}
/**
* @brief Test for integrated passing of the hierarchical basis type
* (needs a lot of boilerplate code to set up manager etc, but its really only the
* `setCombiParametersHierarchicalBasesUniform<T>(params);` and
* `BOOST_TEST(dynamic_cast<T*>(b) != nullptr)` parts that are interesting)
* @tparam T the hierarchical basis class to test (needs to be derived from BasisFunctionBasis,
* serializable, exported (cf. BoostExports.hpp))
*/
template <typename T>
void checkPassingHierarchicalBases(size_t ngroup = 1, size_t nprocs = 1) {
size_t size = ngroup * nprocs + 1;
BOOST_REQUIRE(TestHelper::checkNumMPIProcsAvailable(size));
CommunicatorType comm = TestHelper::getComm(size);
if (comm == MPI_COMM_NULL) {
BOOST_TEST_CHECKPOINT("drop out of test comm");
return;
}
combigrid::Stats::initialize();
theMPISystem()->initWorldReusable(comm, ngroup, nprocs);
DimType dim = 2;
LevelVector lmin(dim, 2);
LevelVector lmax(dim, 5);
WORLD_MANAGER_EXCLUSIVE_SECTION {
ProcessGroupManagerContainer pgroups;
for (int i = 0; i < ngroup; ++i) {
int pgroupRootID(i);
pgroups.emplace_back(std::make_shared<ProcessGroupManager>(pgroupRootID));
}
auto loadmodel = std::unique_ptr<LoadModel>(new LinearLoadModel());
std::vector<BoundaryType> boundary(dim, 2);
CombiMinMaxScheme combischeme(dim, lmin, lmax);
combischeme.createAdaptiveCombischeme();
std::vector<LevelVector> levels = combischeme.getCombiSpaces();
std::vector<combigrid::real> coeffs = combischeme.getCoeffs();
[[maybe_unused]] auto numDOF = printCombiDegreesOfFreedom(levels, boundary);
// create Tasks
TaskContainer tasks;
std::vector<size_t> taskIDs;
for (size_t i = 0; i < levels.size(); i++) {
Task* t = new TaskCount(levels[i], boundary, coeffs[i], loadmodel.get());
tasks.push_back(t);
taskIDs.push_back(t->getID());
}
// create combiparameters
CombiParameters params(dim, lmin, lmax, boundary, levels, coeffs, taskIDs, 2);
params.setParallelization({static_cast<int>(nprocs), 1});
setCombiParametersHierarchicalBasesUniform<T>(params);
// create abstraction for Manager
ProcessManager manager{pgroups, tasks, params, std::move(loadmodel)};
// the combiparameters are sent to all process groups before the
// computations start
manager.updateCombiParameters();
manager.runfirst();
manager.combine();
manager.runnext();
manager.combine();
manager.exit();
BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getGlobalComm()));
}
else {
BOOST_TEST_CHECKPOINT("Worker starts");
ProcessGroupWorker pgroup;
SignalType signal = -1;
// omitting to count RUN_FIRST signal, as it is executed once for every task
while (signal != EXIT) {
BOOST_TEST_CHECKPOINT("Last Successful Worker Signal " + std::to_string(signal));
signal = pgroup.wait();
}
const auto& bases = pgroup.getCombiParameters().getHierarchicalBases();
for (const auto& b : bases) {
BOOST_TEST(dynamic_cast<T*>(b) != nullptr);
}
BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getLocalComm()));
MASTER_EXCLUSIVE_SECTION { BOOST_CHECK(!TestHelper::testStrayMessages(theMPISystem()->getGlobalComm())); }
}
combigrid::Stats::finalize();
MPI_Barrier(comm);
BOOST_CHECK(!TestHelper::testStrayMessages(comm));
}
#ifndef ISGENE // integration tests won't work with ISGENE because of worker magic
#ifndef NDEBUG // in case of a build with asserts, have longer timeout
BOOST_FIXTURE_TEST_SUITE(integration, TestHelper::BarrierAtEnd, *boost::unit_test::timeout(580))
#else
BOOST_FIXTURE_TEST_SUITE(integration, TestHelper::BarrierAtEnd, *boost::unit_test::timeout(480))
#endif // NDEBUG
BOOST_AUTO_TEST_CASE(test_1, *boost::unit_test::tolerance(TestHelper::higherTolerance)) {
auto start = std::chrono::high_resolution_clock::now();
auto rank = TestHelper::getRank(MPI_COMM_WORLD);
for (BoundaryType boundary : std::vector<BoundaryType>({0, 1, 2})) {
for (size_t ngroup : {1, 2, 3, 4}) {
for (size_t nprocs : {1, 2}) {
if (rank == 0)
std::cout << "integration/test_1 " << static_cast<int>(boundary) << " " << ngroup << " "
<< nprocs << std::endl;
BOOST_CHECK_NO_THROW(checkIntegration(ngroup, nprocs, boundary));
MPI_Barrier(MPI_COMM_WORLD);
}
}
for (size_t ngroup : {1, 2}) {
for (size_t nprocs : {3}) {
if (rank == 0)
std::cout << "integration/test_1 " << static_cast<int>(boundary) << " " << ngroup << " "
<< nprocs << std::endl;
BOOST_CHECK_NO_THROW(checkIntegration(ngroup, nprocs, boundary));
MPI_Barrier(MPI_COMM_WORLD);
}
}
for (size_t ngroup : {1, 2}) {
if (boundary > 0) {
for (size_t nprocs : {4}) {
if (rank == 0)
std::cout << "integration/test_1 " << static_cast<int>(boundary) << " " << ngroup << " "
<< nprocs << std::endl;
BOOST_CHECK_NO_THROW(checkIntegration(ngroup, nprocs, boundary));
MPI_Barrier(MPI_COMM_WORLD);
}
}
}
for (size_t ngroup : {1}) {
for (size_t nprocs : {5}) {
if (boundary == 2) {
if (rank == 0)
std::cout << "integration/test_1 " << static_cast<int>(boundary) << " " << ngroup << " "
<< nprocs << std::endl;
BOOST_CHECK_NO_THROW(checkIntegration(ngroup, nprocs, boundary));
MPI_Barrier(MPI_COMM_WORLD);
}
}
}
MPI_Barrier(MPI_COMM_WORLD);
BOOST_CHECK(!TestHelper::testStrayMessages());
}
auto end = std::chrono::high_resolution_clock::now();
auto duration = std::chrono::duration_cast<std::chrono::milliseconds>(end - start);
BOOST_TEST_MESSAGE("time to run all 'integration' tests: " << duration.count()
<< " milliseconds");
}
BOOST_AUTO_TEST_CASE(test_2) { checkPassingHierarchicalBases<HierarchicalHatBasisFunction>(1, 1); }
BOOST_AUTO_TEST_CASE(test_3) { checkPassingHierarchicalBases<FullWeightingBasisFunction>(1, 2); }
BOOST_AUTO_TEST_CASE(test_4) {
checkPassingHierarchicalBases<FullWeightingPeriodicBasisFunction>(2, 2);
}
BOOST_AUTO_TEST_CASE(test_5) { checkPassingHierarchicalBases<BiorthogonalBasisFunction>(1, 4); }
BOOST_AUTO_TEST_CASE(test_6) {
checkPassingHierarchicalBases<BiorthogonalPeriodicBasisFunction>(4, 2);
}
BOOST_AUTO_TEST_CASE(test_7) {
// unit test for downsampleDecomposition
LevelVector lmin{1, 2, 4};
LevelVector lmax{6, 5, 4};
std::vector<IndexVector> decomposition{
{0, 3, 14},
{0},
{0, 1},
};
auto newDecomposition = downsampleDecomposition(decomposition, lmax, lmin, {true, false, false});
BOOST_CHECK(newDecomposition[0] == IndexVector({0, 1, 1}));
BOOST_CHECK(newDecomposition[1] == IndexVector({0}));
BOOST_CHECK(newDecomposition[2] == IndexVector({0, 1}));
}
BOOST_AUTO_TEST_CASE(test_8) {
// unit test for CombiMinMaxSchemeFromFile
LevelVector lmin = {3, 6};
LevelVector lmax = {7, 10};
auto dim = static_cast<DimType>(lmin.size());
std::unique_ptr<CombiMinMaxScheme> scheme(
new CombiMinMaxSchemeFromFile(dim, lmin, lmax, "test_scheme.json"));
auto coeffs = scheme->getCoeffs();
auto levels = scheme->getCombiSpaces();
BOOST_CHECK(std::accumulate(coeffs.begin(), coeffs.end(), 0.) == 1.);
BOOST_CHECK(dynamic_cast<CombiMinMaxSchemeFromFile*>(scheme.get()));
if (auto schemeFromFile = dynamic_cast<CombiMinMaxSchemeFromFile*>(scheme.get())) {
auto pgNumbers = schemeFromFile->getProcessGroupNumbers();
const auto [itMin, itMax] = std::minmax_element(pgNumbers.begin(), pgNumbers.end());
BOOST_CHECK(*itMax == 7);
BOOST_CHECK(*itMin == 0);
auto boundary = std::vector<BoundaryType>(dim, 2);
auto rank = TestHelper::getRank(MPI_COMM_WORLD);
for (size_t taskNo = 0; taskNo < coeffs.size(); ++taskNo) {
BOOST_TEST_CHECKPOINT(std::to_string(rank) + " Last taskNo " + std::to_string(taskNo));
if (pgNumbers[taskNo] == rank) {
auto loadmodel = LinearLoadModel();
TaskCount tc(levels[taskNo], boundary, coeffs[taskNo], &loadmodel);
tc.setID(taskNo);
}
}
}
// clears scheme, may be useful if there are many tasks and
// we only need some allocated on each process group
scheme.reset();
}
BOOST_AUTO_TEST_SUITE_END()
#endif