49 #ifndef Intrepid2_TestUtils_h
50 #define Intrepid2_TestUtils_h
52 #include "Kokkos_Core.hpp"
53 #include "Kokkos_DynRankView.hpp"
61 #include "Teuchos_UnitTestHarness.hpp"
63 static const double TEST_TOLERANCE_TIGHT = 1.e2 * std::numeric_limits<double>::epsilon();
66 template<
typename ScalarType>
67 using ViewType = Kokkos::DynRankView<ScalarType,Kokkos::DefaultExecutionSpace>;
69 template<
typename ScalarType>
70 inline bool valuesAreSmall(ScalarType a, ScalarType b,
double epsilon)
73 return (abs(a) < epsilon) && (abs(b) < epsilon);
76 inline bool approximatelyEqual(
double a,
double b,
double epsilon)
78 const double larger_magnitude = (std::abs(a) < std::abs(b) ? std::abs(b) : std::abs(a));
79 return std::abs(a - b) <= larger_magnitude * epsilon;
82 inline bool essentiallyEqual(
double a,
double b,
double epsilon)
84 const double smaller_magnitude = (std::abs(a) > std::abs(b) ? std::abs(b) : std::abs(a));
85 return std::abs(a - b) <= smaller_magnitude * epsilon;
89 KOKKOS_INLINE_FUNCTION
double fromZeroOne(
double x_zero_one)
91 return x_zero_one * 2.0 - 1.0;
95 KOKKOS_INLINE_FUNCTION
double toZeroOne(
double x_minus_one_one)
97 return (x_minus_one_one + 1.0) / 2.0;
101 KOKKOS_INLINE_FUNCTION
double fromZeroOne_dx(
double dx_zero_one)
103 return dx_zero_one / 2.0;
107 KOKKOS_INLINE_FUNCTION
double toZeroOne_dx(
double dx_minus_one_one)
109 return dx_minus_one_one * 2.0;
112 template<
class DeviceViewType>
113 typename DeviceViewType::HostMirror getHostCopy(
const DeviceViewType &deviceView)
115 typename DeviceViewType::HostMirror hostView = Kokkos::create_mirror(deviceView);
116 Kokkos::deep_copy(hostView, deviceView);
120 template<
class BasisFamily>
121 inline Teuchos::RCP< Intrepid2::Basis<Kokkos::DefaultExecutionSpace,double,double> > getBasisUsingFamily(shards::CellTopology cellTopo, Intrepid2::EFunctionSpace fs,
122 int polyOrder_x,
int polyOrder_y=-1,
int polyOrder_z = -1)
124 using BasisPtr =
typename BasisFamily::BasisPtr;
127 using namespace Intrepid2;
129 if (cellTopo.getBaseKey() == shards::Line<>::key)
131 basis = getLineBasis<BasisFamily>(fs,polyOrder_x);
133 else if (cellTopo.getBaseKey() == shards::Quadrilateral<>::key)
135 INTREPID2_TEST_FOR_EXCEPTION(polyOrder_y < 0, std::invalid_argument,
"polyOrder_y must be specified");
136 basis = getQuadrilateralBasis<BasisFamily>(fs,polyOrder_x,polyOrder_y);
138 else if (cellTopo.getBaseKey() == shards::Triangle<>::key)
140 basis = getTriangleBasis<BasisFamily>(fs,polyOrder_x);
142 else if (cellTopo.getBaseKey() == shards::Hexahedron<>::key)
144 INTREPID2_TEST_FOR_EXCEPTION(polyOrder_y < 0, std::invalid_argument,
"polyOrder_y must be specified");
145 INTREPID2_TEST_FOR_EXCEPTION(polyOrder_z < 0, std::invalid_argument,
"polyOrder_z must be specified");
146 basis = getHexahedronBasis<BasisFamily>(fs,polyOrder_x,polyOrder_y,polyOrder_z);
148 else if (cellTopo.getBaseKey() == shards::Tetrahedron<>::key)
150 basis = getTetrahedronBasis<BasisFamily>(fs, polyOrder_x);
154 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported cell topology");
159 template<
bool defineVertexFunctions>
160 inline Teuchos::RCP< Intrepid2::Basis<Kokkos::DefaultExecutionSpace,double,double> > getHierarchicalBasis(shards::CellTopology cellTopo, Intrepid2::EFunctionSpace fs,
161 int polyOrder_x,
int polyOrder_y=-1,
int polyOrder_z = -1)
163 using ExecSpace = Kokkos::DefaultExecutionSpace;
164 using Scalar = double;
165 using namespace Intrepid2;
174 return getBasisUsingFamily<BasisFamily>(cellTopo, fs, polyOrder_x, polyOrder_y, polyOrder_z);
177 template<
typename ValueType,
class ... DimArgs>
178 inline ViewType<ValueType> getView(
const std::string &label, DimArgs... dims)
180 const bool allocateFadStorage = !std::is_pod<ValueType>::value;
181 constexpr
int maxDerivatives = 3;
182 if (!allocateFadStorage)
184 return ViewType<ValueType>(label,dims...);
188 return ViewType<ValueType>(label,dims...,maxDerivatives+1);
195 template <
typename Po
intValueType>
198 ViewType<PointValueType> inputPoints = getView<PointValueType>(
"line input points",numPoints,1);
199 Kokkos::parallel_for(numPoints, KOKKOS_LAMBDA(
const int i)
201 double x_zero_one = i * 1.0 / (numPoints-1);
202 inputPoints(i,0) = PointValueType(fromZeroOne(x_zero_one));
212 template <
typename Po
intValueType>
215 ViewType<PointValueType> inputPoints = getView<PointValueType>(
"hex input points",numPoints_1D*numPoints_1D*numPoints_1D,3);
216 Kokkos::parallel_for(numPoints_1D, KOKKOS_LAMBDA(
const int i)
218 double x_zero_one = i * 1.0 / (numPoints_1D-1);
219 for (
int j=0; j<numPoints_1D; j++)
221 double y_zero_one = j * 1.0 / (numPoints_1D-1);
222 for (
int k=0; k<numPoints_1D; k++)
224 double z_zero_one = k * 1.0 / (numPoints_1D-1);
225 int pointOrdinal = (i*numPoints_1D+j)*numPoints_1D+k;
226 inputPoints(pointOrdinal,0) = PointValueType(fromZeroOne(x_zero_one));
227 inputPoints(pointOrdinal,1) = PointValueType(fromZeroOne(y_zero_one));
228 inputPoints(pointOrdinal,2) = PointValueType(fromZeroOne(z_zero_one));
240 template <
typename Po
intValueType>
243 ViewType<PointValueType> inputPoints = getView<PointValueType>(
"quad input points",numPoints_1D*numPoints_1D,2);
245 Kokkos::parallel_for(numPoints_1D, KOKKOS_LAMBDA(
const int i)
247 double x_zero_one = i * 1.0 / (numPoints_1D-1);
248 for (
int j=0; j<numPoints_1D; j++)
250 double y_zero_one = j * 1.0 / (numPoints_1D-1);
251 int pointOrdinal = i*numPoints_1D+j;
252 inputPoints(pointOrdinal,0) = PointValueType(fromZeroOne(x_zero_one));
253 inputPoints(pointOrdinal,1) = PointValueType(fromZeroOne(y_zero_one));
264 template <
typename Po
intValueType>
267 const int numPoints = numPointsBase*(numPointsBase+1)*(numPointsBase+2)/6;
268 ViewType<PointValueType> inputPoints = getView<PointValueType>(
"tetrahedron input points",numPoints,3);
269 Kokkos::parallel_for(numPointsBase, KOKKOS_LAMBDA(
const int d0)
275 int pointOrdinalOffset = 0;
276 for (
int i=0; i<d0; i++)
278 for (
int j=0; j<numPointsBase-i; j++)
280 pointOrdinalOffset += (numPointsBase - i - j);
284 double z_zero_one = d0 * 1.0 / (numPointsBase-1);
285 int pointOrdinal = pointOrdinalOffset;
286 for (
int d1=0; d1<numPointsBase-d0; d1++)
288 double y_zero_one = d1 * 1.0 / (numPointsBase-1);
289 for (
int d2=0; d2<numPointsBase-d0-d1; d2++)
291 double x_zero_one = d2 * 1.0 / (numPointsBase-1);
293 inputPoints(pointOrdinal,0) = PointValueType(x_zero_one);
294 inputPoints(pointOrdinal,1) = PointValueType(y_zero_one);
295 inputPoints(pointOrdinal,2) = PointValueType(z_zero_one);
309 template <
typename Po
intValueType>
312 const int numPoints = numPointsBase*(numPointsBase+1)/2;
313 ViewType<PointValueType> inputPoints = getView<PointValueType>(
"triangle input points",numPoints,2);
314 Kokkos::parallel_for(numPointsBase, KOKKOS_LAMBDA(
const int row)
316 int rowPointOrdinalOffset = 0;
317 for (
int i=0; i<row; i++)
319 rowPointOrdinalOffset += (numPointsBase - i);
321 double y_zero_one = row * 1.0 / (numPointsBase-1);
322 for (
int col=0; col<numPointsBase-row; col++)
324 const int pointOrdinal = rowPointOrdinalOffset + col;
325 double x_zero_one = col * 1.0 / (numPointsBase-1);
326 inputPoints(pointOrdinal,0) = PointValueType(x_zero_one);
327 inputPoints(pointOrdinal,1) = PointValueType(y_zero_one);
333 template <
typename Po
intValueType>
334 inline ViewType<PointValueType> getInputPointsView(shards::CellTopology &cellTopo,
int numPoints_1D)
336 if (cellTopo.getBaseKey() == shards::Line<>::key)
338 return lineInputPointsView<PointValueType>(numPoints_1D);
340 else if (cellTopo.getBaseKey() == shards::Quadrilateral<>::key)
342 return quadInputPointsView<PointValueType>(numPoints_1D);
344 else if (cellTopo.getBaseKey() == shards::Hexahedron<>::key)
346 return hexInputPointsView<PointValueType>(numPoints_1D);
348 else if (cellTopo.getBaseKey() == shards::Triangle<>::key)
350 return triInputPointsView<PointValueType>(numPoints_1D);
352 else if (cellTopo.getBaseKey() == shards::Tetrahedron<>::key)
354 return tetInputPointsView<PointValueType>(numPoints_1D);
358 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported Cell Topology");
362 template<
typename OutputValueType>
363 inline ViewType<OutputValueType> getOutputView(Intrepid2::EFunctionSpace fs, Intrepid2::EOperator op,
int basisCardinality,
int numPoints,
int spaceDim)
366 case Intrepid2::FUNCTION_SPACE_HGRAD:
368 case Intrepid2::OPERATOR_VALUE:
369 return getView<OutputValueType>(
"H^1 value output",basisCardinality,numPoints);
370 case Intrepid2::OPERATOR_GRAD:
371 return getView<OutputValueType>(
"H^1 derivative output",basisCardinality,numPoints,spaceDim);
372 case Intrepid2::OPERATOR_D1:
373 case Intrepid2::OPERATOR_D2:
374 case Intrepid2::OPERATOR_D3:
375 case Intrepid2::OPERATOR_D4:
376 case Intrepid2::OPERATOR_D5:
377 case Intrepid2::OPERATOR_D6:
378 case Intrepid2::OPERATOR_D7:
379 case Intrepid2::OPERATOR_D8:
380 case Intrepid2::OPERATOR_D9:
381 case Intrepid2::OPERATOR_D10:
384 return getView<OutputValueType>(
"H^1 derivative output",basisCardinality,numPoints,dkcard);
387 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported op/fs combination");
389 case Intrepid2::FUNCTION_SPACE_HCURL:
391 case Intrepid2::OPERATOR_VALUE:
392 return getView<OutputValueType>(
"H(curl) value output",basisCardinality,numPoints,spaceDim);
393 case Intrepid2::OPERATOR_CURL:
395 return getView<OutputValueType>(
"H(curl) derivative output",basisCardinality,numPoints);
396 else if (spaceDim == 3)
397 return getView<OutputValueType>(
"H(curl) derivative output",basisCardinality,numPoints,spaceDim);
399 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported op/fs combination");
401 case Intrepid2::FUNCTION_SPACE_HDIV:
403 case Intrepid2::OPERATOR_VALUE:
404 return getView<OutputValueType>(
"H(div) value output",basisCardinality,numPoints,spaceDim);
405 case Intrepid2::OPERATOR_DIV:
406 return getView<OutputValueType>(
"H(div) derivative output",basisCardinality,numPoints);
408 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported op/fs combination");
411 case Intrepid2::FUNCTION_SPACE_HVOL:
413 case Intrepid2::OPERATOR_VALUE:
414 return getView<OutputValueType>(
"H(vol) value output",basisCardinality,numPoints);
415 case Intrepid2::OPERATOR_D1:
416 case Intrepid2::OPERATOR_D2:
417 case Intrepid2::OPERATOR_D3:
418 case Intrepid2::OPERATOR_D4:
419 case Intrepid2::OPERATOR_D5:
420 case Intrepid2::OPERATOR_D6:
421 case Intrepid2::OPERATOR_D7:
422 case Intrepid2::OPERATOR_D8:
423 case Intrepid2::OPERATOR_D9:
424 case Intrepid2::OPERATOR_D10:
427 return getView<OutputValueType>(
"H(vol) derivative output",basisCardinality,numPoints,dkcard);
430 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported op/fs combination");
433 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported op/fs combination");
439 inline std::vector< std::vector<int> > getBasisTestCasesUpToDegree(
int spaceDim,
int minDegree,
int polyOrder_x,
int polyOrder_y=-1,
int polyOrder_z = -1)
441 std::vector<int> degrees(spaceDim);
442 degrees[0] = polyOrder_x;
443 if (spaceDim > 1) degrees[1] = polyOrder_y;
444 if (spaceDim > 2) degrees[2] = polyOrder_z;
446 int numCases = degrees[0];
447 for (
unsigned d=1; d<degrees.size(); d++)
449 numCases = numCases * (degrees[d] + 1 - minDegree);
451 std::vector< std::vector<int> > subBasisDegreeTestCases(numCases);
452 for (
int caseOrdinal=0; caseOrdinal<numCases; caseOrdinal++)
454 std::vector<int> subBasisDegrees(degrees.size());
455 int caseRemainder = caseOrdinal;
456 for (
int d=degrees.size()-1; d>=0; d--)
458 int subBasisDegree = caseRemainder % (degrees[d] + 1 - minDegree);
459 caseRemainder = caseRemainder / (degrees[d] + 1 - minDegree);
460 subBasisDegrees[d] = subBasisDegree + minDegree;
462 subBasisDegreeTestCases[caseOrdinal] = subBasisDegrees;
464 return subBasisDegreeTestCases;
467 template <
class ViewType>
468 void testViewFloatingEquality(ViewType &view1, ViewType &view2,
double relTol,
double absTol, Teuchos::FancyOStream &out,
bool &success,
469 std::string view1Name =
"View 1", std::string view2Name =
"View 2")
471 using Scalar =
typename ViewType::value_type;
472 using HostView =
typename ViewType::HostMirror;
475 auto hostView1 = getHostCopy(view1);
476 auto hostView2 = getHostCopy(view2);
479 TEUCHOS_TEST_FOR_EXCEPTION(view1.rank() != view2.rank(), std::invalid_argument,
"views must agree in rank");
480 for (
unsigned i=0; i<view1.rank(); i++)
482 TEUCHOS_TEST_FOR_EXCEPTION(view1.extent_int(i) != view2.extent_int(i), std::invalid_argument,
"views must agree in size in each dimension");
485 if (view1.size() == 0)
return;
487 HostViewIteratorScalar vi1(hostView1);
488 HostViewIteratorScalar vi2(hostView2);
490 double maxDiff = 0.0;
491 double maxRelativeDiff = 0.0;
492 bool differencesFound =
false;
494 bool moreEntries = (vi1.nextIncrementRank() != -1) && (vi2.nextIncrementRank() != -1);
497 Scalar value1 = vi1.get();
498 Scalar value2 = vi2.get();
500 const bool small = valuesAreSmall(value1, value2, absTol);
501 const bool valuesMatch = small || essentiallyEqual(value1, value2, relTol);
505 differencesFound =
true;
507 auto location = vi1.getLocation();
508 out <<
"At location (";
509 for (
unsigned i=0; i<view1.rank(); i++)
512 if (i<view1.rank()-1)
517 out <<
") " << view1Name <<
" value != " << view2Name <<
" value (";
518 out << value1 <<
" != " << value2 <<
"); diff is " << std::abs(value1-value2) << std::endl;
520 maxDiff = std::max(maxDiff, std::abs(value1-value2));
521 double smallerMagnitude = (std::abs(value1) > std::abs(value2) ? std::abs(value2) : std::abs(value1));
522 if (smallerMagnitude > 0)
524 const double relativeDiff = std::abs(value1 - value2) / smallerMagnitude;
525 maxRelativeDiff = std::max(maxRelativeDiff, relativeDiff);
529 moreEntries = (vi1.increment() != -1);
530 moreEntries = moreEntries && (vi2.increment() != -1);
533 if (differencesFound)
535 out <<
"max difference = " << maxDiff << std::endl;
536 out <<
"max relative difference = " << maxRelativeDiff << std::endl;
540 #ifdef HAVE_INTREPID2_SACADO
541 #include <Sacado.hpp>
542 using Sacado_Fad_DFadType = Sacado::Fad::DFad<double>;
543 #define INTREPID2_OUTPUTSCALAR_POINTSCALAR_TEST_INSTANT(GROUP_NAME, TEST_NAME) \
545 TEUCHOS_UNIT_TEST_TEMPLATE_2_INSTANT( GROUP_NAME, TEST_NAME, double, double ) \
547 TEUCHOS_UNIT_TEST_TEMPLATE_2_INSTANT( GROUP_NAME, TEST_NAME, Sacado_Fad_DFadType, Sacado_Fad_DFadType ) \
550 #define INTREPID2_OUTPUTSCALAR_POINTSCALAR_TEST_INSTANT(GROUP_NAME, TEST_NAME) \
552 TEUCHOS_UNIT_TEST_TEMPLATE_2_INSTANT( GROUP_NAME, TEST_NAME, double, double ) \
Header file for the abstract base class Intrepid2::Basis.
KOKKOS_INLINE_FUNCTION ordinal_type getDkCardinality(const EOperator operatorType, const ordinal_type spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.
Stateless class representing a family of basis functions, templated on H(vol) and H(grad) on the line...
Stateless classes that act as factories for two families of hierarchical bases. HierarchicalBasisFami...
Header file to include all Sacado headers that are required if using Intrepid2 with Sacado types.
ViewType< PointValueType > tetInputPointsView(int numPointsBase)
Returns a View containing regularly-spaced points on the tetrahedron.
ViewType< PointValueType > hexInputPointsView(int numPoints_1D)
Returns a View containing equispaced points on the hexahedron.
ViewType< PointValueType > triInputPointsView(int numPointsBase)
Returns a View containing regularly-spaced points on the triangle.
ViewType< PointValueType > quadInputPointsView(int numPoints_1D)
Returns a View containing equispaced points on the quadrilateral.
ViewType< PointValueType > lineInputPointsView(int numPoints)
Returns a View containing equispaced points on the line.
Header function for Intrepid2::Util class and other utility functions.
A family of basis functions, constructed from H(vol) and H(grad) bases on the line.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line.
Basis defining Legendre basis on the line, a polynomial subspace of L^2 (a.k.a. H(vol)) on the line.
A helper class that allows iteration over some part of a Kokkos View, while allowing the calling code...