Skip to content

Commit

Permalink
MueLu: test higher-order line, quad in IntrepidPCoarsenFactory.cpp (#…
Browse files Browse the repository at this point in the history
…13534)

MueLu: test higher-order line, quad in IntrepidPCoarsenFactory.cpp.

We increase the polynomial orders tested in IntrepidPCoarsenFactory.cpp to match the maximal polynomial orders supported by Intrepid2, where this is algorithmically supported: namely, for spectral nodes on the line and quad.  

Note that there is an algorithmic failure in GenerateRepresentativeBasisNodes for equispaced line bases of various orders. For example, for a "low" line basis of order 7, and a "high" line basis of order 10, no "candidate" is generated for the third low-order basis function.  See PR #13534 for some numerical details.
  • Loading branch information
CamelliaDPG authored Oct 17, 2024
1 parent db956c4 commit ffcdd92
Show file tree
Hide file tree
Showing 2 changed files with 68 additions and 31 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -524,7 +524,7 @@ void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer
for (size_t j = 0; j < numFieldsHi; j++)
vmax = std::max(vmax, Teuchos::ScalarTraits<SC>::magnitude(LoValues_host(i, j)));

// 2nd pass: Find all values w/i threshhold of target
// 2nd pass: Find all values w/i threshold of target
for (size_t j = 0; j < numFieldsHi; j++) {
if (Teuchos::ScalarTraits<SC>::magnitude(vmax - LoValues_host(i, j)) < threshold * vmax)
representative_node_candidates[i].push_back(j);
Expand Down
97 changes: 67 additions & 30 deletions packages/muelu/test/unit_tests/IntrepidPCoarsenFactory.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -46,15 +46,21 @@ namespace MueLuTests {

/**** some helper methods and classes by Nate ****/
#ifndef TEST_MORE_COMBINATIONS
static const int MAX_LINE_DEGREE = (Intrepid2::Parameters::MaxOrder < 10) ? Intrepid2::Parameters::MaxOrder : 10;
static const int MAX_QUAD_DEGREE = (Intrepid2::Parameters::MaxOrder < 10) ? Intrepid2::Parameters::MaxOrder : 10;
static const int MAX_HEX_DEGREE = (Intrepid2::Parameters::MaxOrder < 4) ? Intrepid2::Parameters::MaxOrder : 4;
static const int MAX_RANK_COUNT = 4;
static const int MAX_LINE_DEGREE_EQUISPACED = (Intrepid2::Parameters::MaxOrder < 10) ? Intrepid2::Parameters::MaxOrder : 10; // equispaced line basis for p >= 10 leads to failures in GenerateRepresentativeBasisNodes: for p_lo = 7, p_hi = 10, e.g., some low-order basis members have no "candidates" found. (Change this to "11" to see the failure.)
static const int MAX_QUAD_DEGREE_EQUISPACED = (Intrepid2::Parameters::MaxOrder < 10) ? Intrepid2::Parameters::MaxOrder : 10;
static const int MAX_HEX_DEGREE_EQUISPACED = (Intrepid2::Parameters::MaxOrder < 4) ? Intrepid2::Parameters::MaxOrder : 4;
static const int MAX_LINE_DEGREE_SPECTRAL = Intrepid2::Parameters::MaxOrder;
static const int MAX_QUAD_DEGREE_SPECTRAL = Intrepid2::Parameters::MaxOrder;
static const int MAX_HEX_DEGREE_SPECTRAL = (Intrepid2::Parameters::MaxOrder < 4) ? Intrepid2::Parameters::MaxOrder : 4;
static const int MAX_RANK_COUNT = 4;
#else
static const int MAX_LINE_DEGREE = Intrepid2::Parameters::MaxOrder;
static const int MAX_QUAD_DEGREE = Intrepid2::Parameters::MaxOrder;
static const int MAX_HEX_DEGREE = Intrepid2::Parameters::MaxOrder;
static const int MAX_RANK_COUNT = 16;
static const int MAX_LINE_DEGREE_EQUISPACED = 10;
static const int MAX_QUAD_DEGREE_EQUISPACED = 10;
static const int MAX_HEX_DEGREE_EQUISPACED = 10;
static const int MAX_LINE_DEGREE_SPECTRAL = Intrepid2::Parameters::MaxOrder;
static const int MAX_QUAD_DEGREE_SPECTRAL = Intrepid2::Parameters::MaxOrder;
static const int MAX_HEX_DEGREE_SPECTRAL = Intrepid2::Parameters::MaxOrder;
static const int MAX_RANK_COUNT = 16;
#endif

using namespace std;
Expand Down Expand Up @@ -940,51 +946,51 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, FindSeeds_Equispaced_
#include "MueLu_UseShortNames.hpp"
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);
FIND_SEEDS_MACROS;
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE_EQUISPACED, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE_EQUISPACED, Intrepid2::POINTTYPE_EQUISPACED, out, success);
}

TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, FindSeeds_Spectral_Line, Scalar, LocalOrdinal, GlobalOrdinal, Node) {
#include "MueLu_UseShortNames.hpp"
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);
FIND_SEEDS_MACROS;
const Intrepid2::EPointType POINTTYPE_SPECTRAL = static_cast<Intrepid2::EPointType>(1); // Not sure why I have to do this...
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE, POINTTYPE_SPECTRAL, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE, POINTTYPE_SPECTRAL, out, success);
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE_SPECTRAL, POINTTYPE_SPECTRAL, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, LineBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_LINE_DEGREE_SPECTRAL, POINTTYPE_SPECTRAL, out, success);
}

TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, FindSeeds_Equispaced_Quad, Scalar, LocalOrdinal, GlobalOrdinal, Node) {
#include "MueLu_UseShortNames.hpp"
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);
FIND_SEEDS_MACROS;
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE_EQUISPACED, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE_EQUISPACED, Intrepid2::POINTTYPE_EQUISPACED, out, success);
}

TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, FindSeeds_Spectral_Quad, Scalar, LocalOrdinal, GlobalOrdinal, Node) {
#include "MueLu_UseShortNames.hpp"
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);
FIND_SEEDS_MACROS;
const Intrepid2::EPointType POINTTYPE_SPECTRAL = static_cast<Intrepid2::EPointType>(1); // Not sure why I have to do this...
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE, POINTTYPE_SPECTRAL, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE, POINTTYPE_SPECTRAL, out, success);
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE_SPECTRAL, POINTTYPE_SPECTRAL, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, QuadBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE_SPECTRAL, POINTTYPE_SPECTRAL, out, success);
}

TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, FindSeeds_Equispaced_Hex, Scalar, LocalOrdinal, GlobalOrdinal, Node) {
#include "MueLu_UseShortNames.hpp"
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);
FIND_SEEDS_MACROS;
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE_EQUISPACED, Intrepid2::POINTTYPE_EQUISPACED, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_QUAD_DEGREE_EQUISPACED, Intrepid2::POINTTYPE_EQUISPACED, out, success);
}

TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, FindSeeds_Spectral_Hex, Scalar, LocalOrdinal, GlobalOrdinal, Node) {
#include "MueLu_UseShortNames.hpp"
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);
FIND_SEEDS_MACROS;
const Intrepid2::EPointType POINTTYPE_SPECTRAL = static_cast<Intrepid2::EPointType>(1); // Not sure why I have to do this...
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_HEX_DEGREE, POINTTYPE_SPECTRAL, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_HEX_DEGREE, POINTTYPE_SPECTRAL, out, success);
testFindSeedsParallel<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_HEX_DEGREE_SPECTRAL, POINTTYPE_SPECTRAL, out, success);
testFindSeedsSerial<LocalOrdinal, GlobalOrdinal, Node, HexBasis, ES, FC, FCO>(MueLuTests::TestHelpers::Parameters::getLib(), MAX_HEX_DEGREE_SPECTRAL, POINTTYPE_SPECTRAL, out, success);
}

/*********************************************************************************************************************/
Expand Down Expand Up @@ -1717,7 +1723,7 @@ bool test_representative_basis(Teuchos::FancyOStream &out, const std::string &na
for (int d = 0; d < int(lo_DofCoords.extent(1)); d++) {
coords[d] = lo_DofCoords_host(lowOrdinal, d);
}
if (coords[0] == 1.0) {
if (std::fabs(coords[0] - 1.0) < 1e-8) {
int globalOrdinal = loNumbering.getGlobalID(coords);
out << globalOrdinal << ": (";
for (int d = 0; d < int(coords.size()) - 1; d++) {
Expand All @@ -1733,7 +1739,7 @@ bool test_representative_basis(Teuchos::FancyOStream &out, const std::string &na
for (int d = 0; d < int(hi_DofCoords.extent(1)); d++) {
coords[d] = hi_DofCoords_host(highOrdinal, d);
}
if (coords[0] == 1.0) {
if (std::fabs(coords[0] - 1.0) < 1e-8) {
int globalOrdinal = hiNumbering.getGlobalID(coords);
out << globalOrdinal << ": (";
for (int d = 0; d < int(coords.size()) - 1; d++) {
Expand All @@ -1750,9 +1756,20 @@ bool test_representative_basis(Teuchos::FancyOStream &out, const std::string &na

// Correctness Test 1: Make sure that there are no duplicates in the representative lists / no low DOF has no candidates
std::vector<bool> is_candidate(hi_DofCoords.extent(0), false);
bool no_doubles = true;
bool no_doubles = true;
bool no_candidates = false;
{
// DEBUGGING
if ((lowPolyDegree == 7) && (highPolyDegree == 10)) {
cout << "lo = " << lowPolyDegree << "; hi = " << highPolyDegree << std::endl;
}
}

for (int lowOrderDof = 0; no_doubles && lowOrderDof < (int)candidates.size(); lowOrderDof++) {
if (candidates[lowOrderDof].size() == 0) no_doubles = false; // this low DOF has no candidates!
if (candidates[lowOrderDof].size() == 0) {
no_candidates = true; // this low DOF has no candidates!
break;
}
for (int l = 0; l < (int)candidates[lowOrderDof].size(); l++)
if (is_candidate[candidates[lowOrderDof][l]] == false)
is_candidate[candidates[lowOrderDof][l]] = true;
Expand All @@ -1765,6 +1782,10 @@ bool test_representative_basis(Teuchos::FancyOStream &out, const std::string &na
out << "ERROR: " << name << " The 'no duplicates' test fails w/ lo/hi = " << lowPolyDegree << "/" << highPolyDegree << std::endl;
return false;
}
if (no_candidates) {
out << "ERROR: " << name << " The 'no low dof has no candidates' test fails w/ lo/hi = " << lowPolyDegree << "/" << highPolyDegree << std::endl;
return false;
}

// Correctness Test 2: Try 2 elements, in all possible relative orientations, and confirm that the
// "lowest global ordinal" tie-breaker always returns the same thing for both neighbors
Expand Down Expand Up @@ -1958,7 +1979,22 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateRepresentativ
typedef Kokkos::DynRankView<MT, typename Node::device_type> FC;
typedef Intrepid2::Basis_HGRAD_LINE_Cn_FEM<ES, MT, MT> Basis;

bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_LINE_EQUISPACED", Intrepid2::POINTTYPE_EQUISPACED, MAX_LINE_DEGREE);
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_LINE_EQUISPACED", Intrepid2::POINTTYPE_EQUISPACED, MAX_LINE_DEGREE_EQUISPACED);
TEST_EQUALITY(rv, true);
}

TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_LINE_Spectral, Scalar, LocalOrdinal, GlobalOrdinal, Node) {
#include "MueLu_UseShortNames.hpp"
;
MUELU_TESTING_SET_OSTREAM;
MUELU_TESTING_LIMIT_SCOPE(Scalar, GlobalOrdinal, Node);

typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
typedef typename Node::device_type::execution_space ES;
typedef Kokkos::DynRankView<MT, typename Node::device_type> FC;
typedef Intrepid2::Basis_HGRAD_LINE_Cn_FEM<ES, MT, MT> Basis;

bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_LINE_EQUISPACED", Intrepid2::POINTTYPE_WARPBLEND, MAX_LINE_DEGREE_SPECTRAL);
TEST_EQUALITY(rv, true);
}

Expand All @@ -1974,7 +2010,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateRepresentativ
typedef Kokkos::DynRankView<MT, typename Node::device_type> FC;
typedef Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<ES, MT, MT> Basis;

bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_QUAD_EQUISPACED", Intrepid2::POINTTYPE_EQUISPACED, MAX_QUAD_DEGREE);
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_QUAD_EQUISPACED", Intrepid2::POINTTYPE_EQUISPACED, MAX_QUAD_DEGREE_EQUISPACED);
TEST_EQUALITY(rv, true);
}

Expand All @@ -1991,7 +2027,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateRepresentativ
typedef Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<ES, MT, MT> Basis;

const Intrepid2::EPointType POINTTYPE_SPECTRAL = static_cast<Intrepid2::EPointType>(1); // Not sure why I have to do this...
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_QUAD_SPECTRAL", POINTTYPE_SPECTRAL, MAX_QUAD_DEGREE);
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_QUAD_SPECTRAL", POINTTYPE_SPECTRAL, MAX_QUAD_DEGREE_SPECTRAL);
TEST_EQUALITY(rv, true);
}

Expand All @@ -2007,7 +2043,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateRepresentativ
typedef Kokkos::DynRankView<MT, typename Node::device_type> FC;
typedef Intrepid2::Basis_HGRAD_HEX_Cn_FEM<ES, MT, MT> Basis;

bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_HEX_EQUISPACED", Intrepid2::POINTTYPE_EQUISPACED, MAX_HEX_DEGREE);
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_HEX_EQUISPACED", Intrepid2::POINTTYPE_EQUISPACED, MAX_HEX_DEGREE_EQUISPACED);
TEST_EQUALITY(rv, true);
}

Expand All @@ -2024,7 +2060,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateRepresentativ
typedef Intrepid2::Basis_HGRAD_HEX_Cn_FEM<ES, MT, MT> Basis;

const Intrepid2::EPointType POINTTYPE_SPECTRAL = static_cast<Intrepid2::EPointType>(1); // Not sure why I have to do this...
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_HEX_SPECTRAL", POINTTYPE_SPECTRAL, MAX_HEX_DEGREE);
bool rv = test_representative_basis<Scalar, LocalOrdinal, GlobalOrdinal, Node, Basis>(out, " GenerateRepresentativeBasisNodes_HEX_SPECTRAL", POINTTYPE_SPECTRAL, MAX_HEX_DEGREE_SPECTRAL);
TEST_EQUALITY(rv, true);
}

Expand All @@ -2049,7 +2085,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, GenerateLoNodeInHighV
Xpetra::UnderlyingLib lib = MueLuTests::TestHelpers::Parameters::getLib();
RCP<const Teuchos::Comm<int>> comm = TestHelpers::Parameters::getDefaultComm();

int max_degree = MAX_QUAD_DEGREE;
int max_degree = MAX_QUAD_DEGREE_EQUISPACED;
double threshold = 1e-10;
GO gst_invalid = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();

Expand Down Expand Up @@ -2815,6 +2851,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_4_DECL(IntrepidPCoarsenFactory, CreatePreconditioner_
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, CreatePreconditioner_p4, Scalar, LO, GO, Node) \
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_LINE_Equispaced, Scalar, LO, GO, Node) \
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_QUAD_Equispaced, Scalar, LO, GO, Node) \
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_LINE_Spectral, Scalar, LO, GO, Node) \
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_QUAD_Spectral, Scalar, LO, GO, Node) \
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_HEX_Equispaced, Scalar, LO, GO, Node) \
TEUCHOS_UNIT_TEST_TEMPLATE_4_INSTANT(IntrepidPCoarsenFactory, GenerateRepresentativeBasisNodes_HEX_Spectral, Scalar, LO, GO, Node) \
Expand Down

0 comments on commit ffcdd92

Please sign in to comment.