Skip to content

Commit

Permalink
Added small unit-test of HexSCV::determinant.
Browse files Browse the repository at this point in the history
  • Loading branch information
alanw0 committed Oct 14, 2016
1 parent 8b93e18 commit 2571a9c
Show file tree
Hide file tree
Showing 4 changed files with 112 additions and 13 deletions.
19 changes: 6 additions & 13 deletions unit_tests/UnitTest1ElemCoordCheck.C
Original file line number Diff line number Diff line change
Expand Up @@ -9,19 +9,10 @@
#include <stk_mesh/base/FieldBase.hpp>
#include <stk_mesh/base/Field.hpp>
#include <stk_mesh/base/GetEntities.hpp>
#include <stk_io/StkMeshIoBroker.hpp>

void fill_mesh_1_elem_per_proc_hex8(stk::mesh::BulkData& bulk)
{
int nprocs = bulk.parallel_size();
std::string meshSpec = "generated:1x1x"+std::to_string(nprocs);

stk::io::StkMeshIoBroker io(bulk.parallel());
io.set_bulk_data(bulk);
io.add_mesh_database(meshSpec, stk::io::READ_MESH);
io.create_input_mesh();
io.populate_bulk_data();
}
#include "UnitTestUtils.h"

namespace {
unsigned count_locally_owned_elems(const stk::mesh::BulkData& bulk)
{
Expand Down Expand Up @@ -61,6 +52,8 @@ void verify_elems_are_unit_cubes(const stk::mesh::BulkData& bulk)
}
}

}//namespace

TEST(Basic, CheckCoords1Elem)
{
stk::ParallelMachine comm = MPI_COMM_WORLD;
Expand All @@ -69,7 +62,7 @@ TEST(Basic, CheckCoords1Elem)
stk::mesh::MetaData meta(spatialDimension);
stk::mesh::BulkData bulk(meta, comm);

fill_mesh_1_elem_per_proc_hex8(bulk);
unit_test_utils::fill_mesh_1_elem_per_proc_hex8(bulk);

EXPECT_EQ(1u, count_locally_owned_elems(bulk));

Expand Down
76 changes: 76 additions & 0 deletions unit_tests/UnitTestHexSCVDeterminant.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,76 @@
#include <gtest/gtest.h>
#include <limits>

#include <stk_util/parallel/Parallel.hpp>
#include <stk_mesh/base/MetaData.hpp>
#include <stk_mesh/base/BulkData.hpp>
#include <stk_mesh/base/Bucket.hpp>
#include <stk_mesh/base/CoordinateSystems.hpp>
#include <stk_mesh/base/FieldBase.hpp>
#include <stk_mesh/base/Field.hpp>
#include <stk_mesh/base/GetEntities.hpp>

#include <master_element/MasterElement.h>

#include "UnitTestUtils.h"

namespace {

void check_HexSCV_determinant(const stk::mesh::BulkData& bulk)
{
typedef stk::mesh::Field<double,stk::mesh::Cartesian> CoordFieldType;
CoordFieldType* coordField = bulk.mesh_meta_data().get_field<CoordFieldType>(stk::topology::NODE_RANK, "coordinates");
EXPECT_TRUE(coordField != nullptr);

stk::mesh::EntityVector elems;
stk::mesh::get_entities(bulk, stk::topology::ELEM_RANK, elems);

const double tolerance = std::numeric_limits<double>::epsilon();

stk::topology hex8 = stk::topology::HEX_8;

const unsigned spatialDim = bulk.mesh_meta_data().spatial_dimension();
std::vector<double> hex8_node_coords(hex8.num_nodes()*spatialDim, 0.0);
std::vector<double> hex8_scvolumes(hex8.num_nodes(), 0.0);

sierra::nalu::HexSCV hexSCV;
double error[1] = {0};
for(stk::mesh::Entity elem : elems) {
EXPECT_EQ(hex8, bulk.bucket(elem).topology());

const stk::mesh::Entity* nodes = bulk.begin_nodes(elem);
unsigned numNodes = bulk.num_nodes(elem);
unsigned counter = 0;
for(unsigned i=0; i<numNodes; ++i) {
double* nodeCoords = stk::mesh::field_data(*coordField, nodes[i]);

for(unsigned d=0; d<spatialDim; ++d) {
hex8_node_coords[counter++] = nodeCoords[d];
}
}

hexSCV.determinant(1, hex8_node_coords.data(), hex8_scvolumes.data(), error);

EXPECT_EQ(0, error[0]);

for(unsigned i=0; i<hex8.num_nodes(); ++i) {
EXPECT_NEAR(0.125, hex8_scvolumes[i], tolerance);
}
}
}

}//namespace

TEST(HexSCV, determinant)
{
stk::ParallelMachine comm = MPI_COMM_WORLD;

unsigned spatialDimension = 3;
stk::mesh::MetaData meta(spatialDimension);
stk::mesh::BulkData bulk(meta, comm);

unit_test_utils::fill_mesh_1_elem_per_proc_hex8(bulk);

check_HexSCV_determinant(bulk);
}

21 changes: 21 additions & 0 deletions unit_tests/UnitTestUtils.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,21 @@
#include <string>

#include <stk_mesh/base/BulkData.hpp>
#include <stk_io/StkMeshIoBroker.hpp>

namespace unit_test_utils {

void fill_mesh_1_elem_per_proc_hex8(stk::mesh::BulkData& bulk)
{
int nprocs = bulk.parallel_size();
std::string meshSpec = "generated:1x1x"+std::to_string(nprocs);

stk::io::StkMeshIoBroker io(bulk.parallel());
io.set_bulk_data(bulk);
io.add_mesh_database(meshSpec, stk::io::READ_MESH);
io.create_input_mesh();
io.populate_bulk_data();
}

}

9 changes: 9 additions & 0 deletions unit_tests/UnitTestUtils.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,9 @@

#include <stk_mesh/base/BulkData.hpp>

namespace unit_test_utils {

void fill_mesh_1_elem_per_proc_hex8(stk::mesh::BulkData& bulk);

}

0 comments on commit 2571a9c

Please sign in to comment.