Select Git revision
Polynomial1DRootComputer.cpp
test_InterpolateItemArray.cpp 18.60 KiB
#include <catch2/catch_test_macros.hpp>
#include <catch2/matchers/catch_matchers_all.hpp>
#include <language/ast/ASTBuilder.hpp>
#include <language/ast/ASTModulesImporter.hpp>
#include <language/ast/ASTNodeDataTypeBuilder.hpp>
#include <language/ast/ASTNodeExpressionBuilder.hpp>
#include <language/ast/ASTNodeFunctionEvaluationExpressionBuilder.hpp>
#include <language/ast/ASTNodeFunctionExpressionBuilder.hpp>
#include <language/ast/ASTNodeTypeCleaner.hpp>
#include <language/ast/ASTSymbolTableBuilder.hpp>
#include <language/utils/PugsFunctionAdapter.hpp>
#include <language/utils/SymbolTable.hpp>
#include <MeshDataBaseForTests.hpp>
#include <mesh/Connectivity.hpp>
#include <mesh/Mesh.hpp>
#include <mesh/MeshData.hpp>
#include <mesh/MeshDataManager.hpp>
#include <language/utils/InterpolateItemArray.hpp>
#include <pegtl/string_input.hpp>
// clazy:excludeall=non-pod-global-static
TEST_CASE("InterpolateItemArray", "[language]")
{
SECTION("interpolate on all items")
{
auto same_cell_array = [](auto f, auto g) -> bool {
using ItemIdType = typename decltype(f)::index_type;
for (ItemIdType item_id = 0; item_id < f.numberOfItems(); ++item_id) {
for (size_t i = 0; i < f.sizeOfArrays(); ++i) {
if (f[item_id][i] != g[item_id][i]) {
return false;
}
}
}
return true;
};
SECTION("1D")
{
constexpr size_t Dimension = 1;
std::array mesh_list = MeshDataBaseForTests::get().all1DMeshes();
for (const auto& named_mesh : mesh_list) {
SECTION(named_mesh.name())
{
auto mesh_1d = named_mesh.mesh();
auto xj = MeshDataManager::instance().getMeshData(*mesh_1d).xj();
std::string_view data = R"(
import math;
let scalar_affine_1d: R^1 -> R, x -> 2*x[0] + 2;
let scalar_non_linear_1d: R^1 -> R, x -> 2 * exp(x[0]) + 3;
)";
TAO_PEGTL_NAMESPACE::string_input input{data, "test.pgs"};
auto ast = ASTBuilder::build(input);
ASTModulesImporter{*ast};
ASTNodeTypeCleaner<language::import_instruction>{*ast};
ASTSymbolTableBuilder{*ast};
ASTNodeDataTypeBuilder{*ast};
ASTNodeTypeCleaner<language::var_declaration>{*ast};
ASTNodeTypeCleaner<language::fct_declaration>{*ast};
ASTNodeExpressionBuilder{*ast};
std::shared_ptr<SymbolTable> symbol_table = ast->m_symbol_table;
TAO_PEGTL_NAMESPACE::position position{TAO_PEGTL_NAMESPACE::internal::iterator{"fixture"}, "fixture"};
position.byte = data.size(); // ensure that variables are declared at this point
std::vector<FunctionSymbolId> function_symbol_id_list;
{
auto [i_symbol, found] = symbol_table->find("scalar_affine_1d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
{
auto [i_symbol, found] = symbol_table->find("scalar_non_linear_1d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
CellArray<double> cell_array{mesh_1d->connectivity(), 2};
parallel_for(
cell_array.numberOfItems(), PUGS_LAMBDA(const CellId cell_id) {
const TinyVector<Dimension>& x = xj[cell_id];
cell_array[cell_id][0] = 2 * x[0] + 2;
cell_array[cell_id][1] = 2 * exp(x[0]) + 3;
});
CellArray<const double> interpolate_array =
InterpolateItemArray<double(TinyVector<Dimension>)>::interpolate(function_symbol_id_list, xj);
REQUIRE(same_cell_array(cell_array, interpolate_array));
}
}
}
SECTION("2D")
{
constexpr size_t Dimension = 2;
std::array mesh_list = MeshDataBaseForTests::get().all2DMeshes();
for (const auto& named_mesh : mesh_list) {
SECTION(named_mesh.name())
{
auto mesh_2d = named_mesh.mesh();
auto xj = MeshDataManager::instance().getMeshData(*mesh_2d).xj();
std::string_view data = R"(
import math;
let scalar_affine_2d: R^2 -> R, x -> 2*x[0] + 3*x[1] + 2;
let scalar_non_linear_2d: R^2 -> R, x -> 2*exp(x[0])*sin(x[1])+3;
)";
TAO_PEGTL_NAMESPACE::string_input input{data, "test.pgs"};
auto ast = ASTBuilder::build(input);
ASTModulesImporter{*ast};
ASTNodeTypeCleaner<language::import_instruction>{*ast};
ASTSymbolTableBuilder{*ast};
ASTNodeDataTypeBuilder{*ast};
ASTNodeTypeCleaner<language::var_declaration>{*ast};
ASTNodeTypeCleaner<language::fct_declaration>{*ast};
ASTNodeExpressionBuilder{*ast};
std::shared_ptr<SymbolTable> symbol_table = ast->m_symbol_table;
TAO_PEGTL_NAMESPACE::position position{TAO_PEGTL_NAMESPACE::internal::iterator{"fixture"}, "fixture"};
position.byte = data.size(); // ensure that variables are declared at this point
std::vector<FunctionSymbolId> function_symbol_id_list;
{
auto [i_symbol, found] = symbol_table->find("scalar_affine_2d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
{
auto [i_symbol, found] = symbol_table->find("scalar_non_linear_2d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
CellArray<double> cell_array{mesh_2d->connectivity(), 2};
parallel_for(
cell_array.numberOfItems(), PUGS_LAMBDA(const CellId cell_id) {
const TinyVector<Dimension>& x = xj[cell_id];
cell_array[cell_id][0] = 2 * x[0] + 3 * x[1] + 2;
cell_array[cell_id][1] = 2 * exp(x[0]) * sin(x[1]) + 3;
});
CellArray<const double> interpolate_array =
InterpolateItemArray<double(TinyVector<Dimension>)>::interpolate(function_symbol_id_list, xj);
REQUIRE(same_cell_array(cell_array, interpolate_array));
}
}
}
SECTION("3D")
{
constexpr size_t Dimension = 3;
std::array mesh_list = MeshDataBaseForTests::get().all3DMeshes();
for (const auto& named_mesh : mesh_list) {
SECTION(named_mesh.name())
{
auto mesh_3d = named_mesh.mesh();
auto xj = MeshDataManager::instance().getMeshData(*mesh_3d).xj();
std::string_view data = R"(
import math;
let scalar_affine_3d: R^3 -> R, x -> 2 * x[0] + 3 * x[1] + 2 * x[2] - 1;
let scalar_non_linear_3d: R^3 -> R, x -> 2 * exp(x[0]) * sin(x[1]) * x[2] + 3;
)";
TAO_PEGTL_NAMESPACE::string_input input{data, "test.pgs"};
auto ast = ASTBuilder::build(input);
ASTModulesImporter{*ast};
ASTNodeTypeCleaner<language::import_instruction>{*ast};
ASTSymbolTableBuilder{*ast};
ASTNodeDataTypeBuilder{*ast};
ASTNodeTypeCleaner<language::var_declaration>{*ast};
ASTNodeTypeCleaner<language::fct_declaration>{*ast};
ASTNodeExpressionBuilder{*ast};
std::shared_ptr<SymbolTable> symbol_table = ast->m_symbol_table;
TAO_PEGTL_NAMESPACE::position position{TAO_PEGTL_NAMESPACE::internal::iterator{"fixture"}, "fixture"};
position.byte = data.size(); // ensure that variables are declared at this point
std::vector<FunctionSymbolId> function_symbol_id_list;
{
auto [i_symbol, found] = symbol_table->find("scalar_affine_3d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
{
auto [i_symbol, found] = symbol_table->find("scalar_non_linear_3d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
CellArray<double> cell_array{mesh_3d->connectivity(), 2};
parallel_for(
cell_array.numberOfItems(), PUGS_LAMBDA(const CellId cell_id) {
const TinyVector<Dimension>& x = xj[cell_id];
cell_array[cell_id][0] = 2 * x[0] + 3 * x[1] + 2 * x[2] - 1;
cell_array[cell_id][1] = 2 * exp(x[0]) * sin(x[1]) * x[2] + 3;
});
CellArray<const double> interpolate_array =
InterpolateItemArray<double(TinyVector<Dimension>)>::interpolate(function_symbol_id_list, xj);
REQUIRE(same_cell_array(cell_array, interpolate_array));
}
}
}
}
SECTION("interpolate on items list")
{
auto same_cell_value = [](auto interpolated, auto reference) -> bool {
for (size_t i = 0; i < interpolated.numberOfRows(); ++i) {
for (size_t j = 0; j < interpolated.numberOfColumns(); ++j) {
if (interpolated[i][j] != reference[i][j]) {
return false;
}
}
}
return true;
};
SECTION("1D")
{
constexpr size_t Dimension = 1;
std::array mesh_list = MeshDataBaseForTests::get().all1DMeshes();
for (const auto& named_mesh : mesh_list) {
SECTION(named_mesh.name())
{
auto mesh_1d = named_mesh.mesh();
auto xj = MeshDataManager::instance().getMeshData(*mesh_1d).xj();
Array<const CellId> cell_id_list = [&] {
Array<CellId> cell_ids{mesh_1d->numberOfCells() / 2};
for (size_t i_cell = 0; i_cell < cell_ids.size(); ++i_cell) {
cell_ids[i_cell] = static_cast<CellId>(2 * i_cell);
}
return cell_ids;
}();
std::string_view data = R"(
import math;
let scalar_affine_1d: R^1 -> R, x -> 2*x[0] + 2;
let scalar_non_linear_1d: R^1 -> R, x -> 2 * exp(x[0]) + 3;
)";
TAO_PEGTL_NAMESPACE::string_input input{data, "test.pgs"};
auto ast = ASTBuilder::build(input);
ASTModulesImporter{*ast};
ASTNodeTypeCleaner<language::import_instruction>{*ast};
ASTSymbolTableBuilder{*ast};
ASTNodeDataTypeBuilder{*ast};
ASTNodeTypeCleaner<language::var_declaration>{*ast};
ASTNodeTypeCleaner<language::fct_declaration>{*ast};
ASTNodeExpressionBuilder{*ast};
std::shared_ptr<SymbolTable> symbol_table = ast->m_symbol_table;
TAO_PEGTL_NAMESPACE::position position{TAO_PEGTL_NAMESPACE::internal::iterator{"fixture"}, "fixture"};
position.byte = data.size(); // ensure that variables are declared at this point
std::vector<FunctionSymbolId> function_symbol_id_list;
{
auto [i_symbol, found] = symbol_table->find("scalar_affine_1d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
{
auto [i_symbol, found] = symbol_table->find("scalar_non_linear_1d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
Table<double> cell_array{cell_id_list.size(), 2};
parallel_for(
cell_id_list.size(), PUGS_LAMBDA(const size_t i) {
const TinyVector<Dimension>& x = xj[cell_id_list[i]];
cell_array[i][0] = 2 * x[0] + 2;
cell_array[i][1] = 2 * exp(x[0]) + 3;
});
Table<const double> interpolate_array =
InterpolateItemArray<double(TinyVector<Dimension>)>::interpolate(function_symbol_id_list, xj, cell_id_list);
REQUIRE(same_cell_value(cell_array, interpolate_array));
}
}
}
SECTION("2D")
{
constexpr size_t Dimension = 2;
std::array mesh_list = MeshDataBaseForTests::get().all2DMeshes();
for (const auto& named_mesh : mesh_list) {
SECTION(named_mesh.name())
{
auto mesh_2d = named_mesh.mesh();
auto xj = MeshDataManager::instance().getMeshData(*mesh_2d).xj();
Array<CellId> cell_id_list{mesh_2d->numberOfCells() / 2};
for (size_t i_cell = 0; i_cell < cell_id_list.size(); ++i_cell) {
cell_id_list[i_cell] = static_cast<CellId>(2 * i_cell);
}
std::string_view data = R"(
import math;
let scalar_affine_2d: R^2 -> R, x -> 2*x[0] + 3*x[1] + 2;
let scalar_non_linear_2d: R^2 -> R, x -> 2*exp(x[0])*sin(x[1])+3;
)";
TAO_PEGTL_NAMESPACE::string_input input{data, "test.pgs"};
auto ast = ASTBuilder::build(input);
ASTModulesImporter{*ast};
ASTNodeTypeCleaner<language::import_instruction>{*ast};
ASTSymbolTableBuilder{*ast};
ASTNodeDataTypeBuilder{*ast};
ASTNodeTypeCleaner<language::var_declaration>{*ast};
ASTNodeTypeCleaner<language::fct_declaration>{*ast};
ASTNodeExpressionBuilder{*ast};
std::shared_ptr<SymbolTable> symbol_table = ast->m_symbol_table;
TAO_PEGTL_NAMESPACE::position position{TAO_PEGTL_NAMESPACE::internal::iterator{"fixture"}, "fixture"};
position.byte = data.size(); // ensure that variables are declared at this point
std::vector<FunctionSymbolId> function_symbol_id_list;
{
auto [i_symbol, found] = symbol_table->find("scalar_affine_2d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
{
auto [i_symbol, found] = symbol_table->find("scalar_non_linear_2d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
Table<double> cell_array{cell_id_list.size(), 2};
parallel_for(
cell_id_list.size(), PUGS_LAMBDA(const size_t i) {
const TinyVector<Dimension>& x = xj[cell_id_list[i]];
cell_array[i][0] = 2 * x[0] + 3 * x[1] + 2;
cell_array[i][1] = 2 * exp(x[0]) * sin(x[1]) + 3;
});
Table<const double> interpolate_array =
InterpolateItemArray<double(TinyVector<Dimension>)>::interpolate(function_symbol_id_list, xj, cell_id_list);
REQUIRE(same_cell_value(cell_array, interpolate_array));
}
}
}
SECTION("3D")
{
constexpr size_t Dimension = 3;
std::array mesh_list = MeshDataBaseForTests::get().all3DMeshes();
for (const auto& named_mesh : mesh_list) {
SECTION(named_mesh.name())
{
auto mesh_3d = named_mesh.mesh();
auto xj = MeshDataManager::instance().getMeshData(*mesh_3d).xj();
Array<CellId> cell_id_list{mesh_3d->numberOfCells() / 2};
for (size_t i_cell = 0; i_cell < cell_id_list.size(); ++i_cell) {
cell_id_list[i_cell] = static_cast<CellId>(2 * i_cell);
}
std::string_view data = R"(
import math;
let scalar_affine_3d: R^3 -> R, x -> 2 * x[0] + 3 * x[1] + 2 * x[2] - 1;
let scalar_non_linear_3d: R^3 -> R, x -> 2 * exp(x[0]) * sin(x[1]) * x[2] + 3;
)";
TAO_PEGTL_NAMESPACE::string_input input{data, "test.pgs"};
auto ast = ASTBuilder::build(input);
ASTModulesImporter{*ast};
ASTNodeTypeCleaner<language::import_instruction>{*ast};
ASTSymbolTableBuilder{*ast};
ASTNodeDataTypeBuilder{*ast};
ASTNodeTypeCleaner<language::var_declaration>{*ast};
ASTNodeTypeCleaner<language::fct_declaration>{*ast};
ASTNodeExpressionBuilder{*ast};
std::shared_ptr<SymbolTable> symbol_table = ast->m_symbol_table;
TAO_PEGTL_NAMESPACE::position position{TAO_PEGTL_NAMESPACE::internal::iterator{"fixture"}, "fixture"};
position.byte = data.size(); // ensure that variables are declared at this point
std::vector<FunctionSymbolId> function_symbol_id_list;
{
auto [i_symbol, found] = symbol_table->find("scalar_affine_3d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
{
auto [i_symbol, found] = symbol_table->find("scalar_non_linear_3d", position);
REQUIRE(found);
REQUIRE(i_symbol->attributes().dataType() == ASTNodeDataType::function_t);
function_symbol_id_list.push_back(
FunctionSymbolId(std::get<uint64_t>(i_symbol->attributes().value()), symbol_table));
}
Table<double> cell_array{cell_id_list.size(), 2};
parallel_for(
cell_id_list.size(), PUGS_LAMBDA(const size_t i) {
const TinyVector<Dimension>& x = xj[cell_id_list[i]];
cell_array[i][0] = 2 * x[0] + 3 * x[1] + 2 * x[2] - 1;
cell_array[i][1] = 2 * exp(x[0]) * sin(x[1]) * x[2] + 3;
});
Table<const double> interpolate_array =
InterpolateItemArray<double(TinyVector<Dimension>)>::interpolate(function_symbol_id_list, xj, cell_id_list);
REQUIRE(same_cell_value(cell_array, interpolate_array));
}
}
}
}
}