Skip to content

Commit d54ac14

Browse files
authored
Merge pull request #17 from alanw0/master
Basic unit-test infrastructure.
2 parents 636a674 + 2571a9c commit d54ac14

File tree

6 files changed

+219
-0
lines changed

6 files changed

+219
-0
lines changed

CMakeLists.txt

Lines changed: 11 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -98,6 +98,8 @@ add_library (nalu ${SOURCE} ${HEADER})
9898
target_link_libraries(nalu ${Trilinos_LIBRARIES})
9999
target_link_libraries(nalu ${YAML_LIBRARY})
100100

101+
file (GLOB UNIT_TESTS_SOURCES unit_tests/*.C)
102+
101103
set(nalu_ex_name "naluX")
102104
message("CMAKE_BUILD_TYPE = ${CMAKE_BUILD_TYPE}")
103105
if (CMAKE_BUILD_TYPE STREQUAL "DEBUG")
@@ -108,4 +110,13 @@ endif()
108110

109111
add_executable(${nalu_ex_name} nalu.C)
110112
target_link_libraries(${nalu_ex_name} nalu)
113+
114+
set(utest_ex_name "unittestX")
115+
if (CMAKE_BUILD_TYPE STREQUAL "DEBUG")
116+
set(utest_ex_name "unittestXd")
117+
endif()
118+
119+
add_executable(${utest_ex_name} unit_tests.C ${UNIT_TESTS_SOURCES})
120+
target_link_libraries(${utest_ex_name} nalu)
121+
111122
MESSAGE("\nAnd CMake says...:")

unit_tests.C

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,31 @@
1+
/*------------------------------------------------------------------------*/
2+
/* Copyright 2014 Sandia Corporation. */
3+
/* This software is released under the license detailed */
4+
/* in the file, LICENSE, which is located in the top-level Nalu */
5+
/* directory structure */
6+
/*------------------------------------------------------------------------*/
7+
8+
#include <gtest/gtest.h> // for InitGoogleTest, etc
9+
#include <mpi.h> // for MPI_Comm_rank, MPI_Finalize, etc
10+
11+
#include "include/NaluEnv.h"
12+
13+
int gl_argc = 0;
14+
char** gl_argv = 0;
15+
16+
int main(int argc, char **argv)
17+
{
18+
MPI_Init(&argc, &argv);
19+
//NaluEnv will call MPI_Finalize for us.
20+
sierra::nalu::NaluEnv::self();
21+
22+
testing::InitGoogleTest(&argc, argv);
23+
24+
gl_argc = argc;
25+
gl_argv = argv;
26+
27+
int returnVal = RUN_ALL_TESTS();
28+
29+
return returnVal;
30+
}
31+
Lines changed: 71 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,71 @@
1+
#include <gtest/gtest.h>
2+
#include <limits>
3+
4+
#include <stk_util/parallel/Parallel.hpp>
5+
#include <stk_mesh/base/MetaData.hpp>
6+
#include <stk_mesh/base/BulkData.hpp>
7+
#include <stk_mesh/base/Bucket.hpp>
8+
#include <stk_mesh/base/CoordinateSystems.hpp>
9+
#include <stk_mesh/base/FieldBase.hpp>
10+
#include <stk_mesh/base/Field.hpp>
11+
#include <stk_mesh/base/GetEntities.hpp>
12+
13+
#include "UnitTestUtils.h"
14+
15+
namespace {
16+
17+
unsigned count_locally_owned_elems(const stk::mesh::BulkData& bulk)
18+
{
19+
return stk::mesh::count_selected_entities(bulk.mesh_meta_data().locally_owned_part(),
20+
bulk.buckets(stk::topology::ELEM_RANK));
21+
}
22+
23+
void verify_elems_are_unit_cubes(const stk::mesh::BulkData& bulk)
24+
{
25+
typedef stk::mesh::Field<double,stk::mesh::Cartesian> CoordFieldType;
26+
CoordFieldType* coordField = bulk.mesh_meta_data().get_field<CoordFieldType>(stk::topology::NODE_RANK, "coordinates");
27+
EXPECT_TRUE(coordField != nullptr);
28+
29+
stk::mesh::EntityVector elems;
30+
stk::mesh::get_entities(bulk, stk::topology::ELEM_RANK, elems);
31+
32+
const double tolerance = std::numeric_limits<double>::epsilon();
33+
34+
for(stk::mesh::Entity elem : elems) {
35+
double minX = 999.99, minY = 999.99, minZ = 999.99;
36+
double maxX = 0, maxY = 0, maxZ = 0;
37+
const stk::mesh::Entity* nodes = bulk.begin_nodes(elem);
38+
unsigned numNodes = bulk.num_nodes(elem);
39+
for(unsigned i=0; i<numNodes; ++i) {
40+
double* nodeCoords = stk::mesh::field_data(*coordField, nodes[i]);
41+
minX = std::min(minX, nodeCoords[0]);
42+
minY = std::min(minY, nodeCoords[1]);
43+
minZ = std::min(minZ, nodeCoords[2]);
44+
maxX = std::max(maxX, nodeCoords[0]);
45+
maxY = std::max(maxY, nodeCoords[1]);
46+
maxZ = std::max(maxZ, nodeCoords[2]);
47+
}
48+
49+
EXPECT_NEAR(1.0, (maxX-minX), tolerance);
50+
EXPECT_NEAR(1.0, (maxY-minY), tolerance);
51+
EXPECT_NEAR(1.0, (maxZ-minZ), tolerance);
52+
}
53+
}
54+
55+
}//namespace
56+
57+
TEST(Basic, CheckCoords1Elem)
58+
{
59+
stk::ParallelMachine comm = MPI_COMM_WORLD;
60+
61+
unsigned spatialDimension = 3;
62+
stk::mesh::MetaData meta(spatialDimension);
63+
stk::mesh::BulkData bulk(meta, comm);
64+
65+
unit_test_utils::fill_mesh_1_elem_per_proc_hex8(bulk);
66+
67+
EXPECT_EQ(1u, count_locally_owned_elems(bulk));
68+
69+
verify_elems_are_unit_cubes(bulk);
70+
}
71+
Lines changed: 76 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,76 @@
1+
#include <gtest/gtest.h>
2+
#include <limits>
3+
4+
#include <stk_util/parallel/Parallel.hpp>
5+
#include <stk_mesh/base/MetaData.hpp>
6+
#include <stk_mesh/base/BulkData.hpp>
7+
#include <stk_mesh/base/Bucket.hpp>
8+
#include <stk_mesh/base/CoordinateSystems.hpp>
9+
#include <stk_mesh/base/FieldBase.hpp>
10+
#include <stk_mesh/base/Field.hpp>
11+
#include <stk_mesh/base/GetEntities.hpp>
12+
13+
#include <master_element/MasterElement.h>
14+
15+
#include "UnitTestUtils.h"
16+
17+
namespace {
18+
19+
void check_HexSCV_determinant(const stk::mesh::BulkData& bulk)
20+
{
21+
typedef stk::mesh::Field<double,stk::mesh::Cartesian> CoordFieldType;
22+
CoordFieldType* coordField = bulk.mesh_meta_data().get_field<CoordFieldType>(stk::topology::NODE_RANK, "coordinates");
23+
EXPECT_TRUE(coordField != nullptr);
24+
25+
stk::mesh::EntityVector elems;
26+
stk::mesh::get_entities(bulk, stk::topology::ELEM_RANK, elems);
27+
28+
const double tolerance = std::numeric_limits<double>::epsilon();
29+
30+
stk::topology hex8 = stk::topology::HEX_8;
31+
32+
const unsigned spatialDim = bulk.mesh_meta_data().spatial_dimension();
33+
std::vector<double> hex8_node_coords(hex8.num_nodes()*spatialDim, 0.0);
34+
std::vector<double> hex8_scvolumes(hex8.num_nodes(), 0.0);
35+
36+
sierra::nalu::HexSCV hexSCV;
37+
double error[1] = {0};
38+
for(stk::mesh::Entity elem : elems) {
39+
EXPECT_EQ(hex8, bulk.bucket(elem).topology());
40+
41+
const stk::mesh::Entity* nodes = bulk.begin_nodes(elem);
42+
unsigned numNodes = bulk.num_nodes(elem);
43+
unsigned counter = 0;
44+
for(unsigned i=0; i<numNodes; ++i) {
45+
double* nodeCoords = stk::mesh::field_data(*coordField, nodes[i]);
46+
47+
for(unsigned d=0; d<spatialDim; ++d) {
48+
hex8_node_coords[counter++] = nodeCoords[d];
49+
}
50+
}
51+
52+
hexSCV.determinant(1, hex8_node_coords.data(), hex8_scvolumes.data(), error);
53+
54+
EXPECT_EQ(0, error[0]);
55+
56+
for(unsigned i=0; i<hex8.num_nodes(); ++i) {
57+
EXPECT_NEAR(0.125, hex8_scvolumes[i], tolerance);
58+
}
59+
}
60+
}
61+
62+
}//namespace
63+
64+
TEST(HexSCV, determinant)
65+
{
66+
stk::ParallelMachine comm = MPI_COMM_WORLD;
67+
68+
unsigned spatialDimension = 3;
69+
stk::mesh::MetaData meta(spatialDimension);
70+
stk::mesh::BulkData bulk(meta, comm);
71+
72+
unit_test_utils::fill_mesh_1_elem_per_proc_hex8(bulk);
73+
74+
check_HexSCV_determinant(bulk);
75+
}
76+

unit_tests/UnitTestUtils.C

Lines changed: 21 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,21 @@
1+
#include <string>
2+
3+
#include <stk_mesh/base/BulkData.hpp>
4+
#include <stk_io/StkMeshIoBroker.hpp>
5+
6+
namespace unit_test_utils {
7+
8+
void fill_mesh_1_elem_per_proc_hex8(stk::mesh::BulkData& bulk)
9+
{
10+
int nprocs = bulk.parallel_size();
11+
std::string meshSpec = "generated:1x1x"+std::to_string(nprocs);
12+
13+
stk::io::StkMeshIoBroker io(bulk.parallel());
14+
io.set_bulk_data(bulk);
15+
io.add_mesh_database(meshSpec, stk::io::READ_MESH);
16+
io.create_input_mesh();
17+
io.populate_bulk_data();
18+
}
19+
20+
}
21+

unit_tests/UnitTestUtils.h

Lines changed: 9 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,9 @@
1+
2+
#include <stk_mesh/base/BulkData.hpp>
3+
4+
namespace unit_test_utils {
5+
6+
void fill_mesh_1_elem_per_proc_hex8(stk::mesh::BulkData& bulk);
7+
8+
}
9+

0 commit comments

Comments
 (0)