Thanks to visit codestin.com
Credit goes to github.com

Skip to content

Commit 97f633c

Browse files
committed
Simplify logic in fuzzy_equals.h
Developers should overload `libMesh::l1_norm` and `libMesh::l1_norm_diff` for their numeric types and then they automatically get correct application of `libMesh::*fuzzy_equals`
1 parent 65d7c4b commit 97f633c

File tree

3 files changed

+68
-110
lines changed

3 files changed

+68
-110
lines changed

include/numerics/type_vector.h

Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1186,6 +1186,21 @@ outer_product(const TypeVector<T> & a, const T2 & b)
11861186

11871187
return ret;
11881188
}
1189+
1190+
template <typename T>
1191+
auto
1192+
l1_norm(const TypeVector<T> & var) -> decltype(var.l1_norm())
1193+
{
1194+
return var.l1_norm();
1195+
}
1196+
1197+
template <typename T>
1198+
auto
1199+
l1_norm_diff(const TypeVector<T> & var) -> decltype(var.l1_norm_diff())
1200+
{
1201+
return var.l1_norm_diff();
1202+
}
1203+
11891204
} // namespace libMesh
11901205

11911206
namespace std

include/utils/fuzzy_equals.h

Lines changed: 50 additions & 110 deletions
Original file line numberDiff line numberDiff line change
@@ -23,57 +23,61 @@
2323
#include "libmesh/compare_types.h"
2424
#include "libmesh/int_range.h"
2525

26+
#ifdef LIBMESH_HAVE_METAPHYSICL
27+
#include "metaphysicl/raw_type.h"
28+
#endif
29+
2630
namespace libMesh
2731
{
2832
/**
29-
* Function to check whether two variables are equal within an absolute tolerance
30-
* @param var1 The first variable to be checked
31-
* @param var2 The second variable to be checked
32-
* @param tol The tolerance to be used
33-
* @return true if var1 and var2 are equal within tol
33+
* Computes the L1 norm
3434
*/
35-
template <typename T,
36-
typename T2,
37-
typename std::enable_if<ScalarTraits<T>::value && ScalarTraits<T2>::value, int>::type = 0>
38-
bool absolute_fuzzy_equals(const T & var1, const T2 & var2, const Real tol = TOLERANCE * TOLERANCE);
35+
template <typename T>
36+
auto
37+
l1_norm(const T & var) -> decltype(std::abs(T{}))
38+
{
39+
return std::abs(var);
40+
}
41+
42+
template <typename T>
43+
auto
44+
l1_norm(const std::complex<T> & var) -> decltype(std::abs(T{}))
45+
{
46+
return std::abs(var.real()) + std::abs(var.imag());
47+
}
3948

4049
/**
41-
* Function to check whether both the real and imaginary parts of two complex variables are equal
42-
* within an absolute tolerance
43-
* @param var1 The first complex variable to be checked
44-
* @param var2 The second complex variable to be checked
45-
* @param tol The real tolerance to be used for comparing both real and imaginary parts
46-
* @return true if both real and imginary parts are equal within tol
50+
* Computes the L1 norm of the diff between \p var1 and \p var2
4751
*/
4852
template <typename T, typename T2>
49-
bool
50-
absolute_fuzzy_equals(const std::complex<T> & var1,
51-
const std::complex<T2> & var2,
52-
const Real tol = TOLERANCE * TOLERANCE)
53+
auto
54+
l1_norm_diff(const T & var1,
55+
const T2 & var2) -> decltype(std::abs(typename CompareTypes<T, T2>::supertype{}))
5356
{
54-
return absolute_fuzzy_equals(var1.real(), var2.real(), tol) &&
55-
absolute_fuzzy_equals(var1.imag(), var2.imag(), tol);
57+
return std::abs(var1 - var2);
5658
}
5759

58-
#ifdef LIBMESH_HAVE_METAPHYSICL
59-
60-
#include "metaphysicl/raw_type.h"
60+
template <typename T, typename T2>
61+
auto
62+
l1_norm_diff(const std::complex<T> & var1, const std::complex<T2> & var2)
63+
-> decltype(std::abs(typename CompareTypes<T, T2>::supertype{}))
64+
{
65+
return std::abs(var1.real() - var2.real()) + std::abs(var1.imag() - var2.imag());
66+
}
6167

68+
#ifdef LIBMESH_HAVE_METAPHYSICL
6269
/**
6370
* Function to check whether two variables are equal within an absolute tolerance
6471
* @param var1 The first variable to be checked
6572
* @param var2 The second variable to be checked
6673
* @param tol The tolerance to be used
6774
* @return true if var1 and var2 are equal within tol
6875
*/
69-
template <typename T,
70-
typename T2,
71-
typename std::enable_if<ScalarTraits<T>::value && ScalarTraits<T2>::value, int>::type>
76+
template <typename T, typename T2>
7277
bool
73-
absolute_fuzzy_equals(const T & var1, const T2 & var2, const Real tol)
78+
absolute_fuzzy_equals(const T & var1, const T2 & var2, const Real tol = TOLERANCE * TOLERANCE)
7479
{
75-
return (std::abs(MetaPhysicL::raw_value(var1) - MetaPhysicL::raw_value(var2)) <=
76-
MetaPhysicL::raw_value(tol));
80+
return MetaPhysicL::raw_value(l1_norm_diff(var1, var2)) <= tol;
7781
}
7882

7983
/**
@@ -87,54 +91,24 @@ template <typename T, typename T2>
8791
bool
8892
relative_fuzzy_equals(const T & var1, const T2 & var2, const Real tol = TOLERANCE * TOLERANCE)
8993
{
90-
if constexpr (ScalarTraits<T>::value || TensorTools::MathWrapperTraits<T>::value)
91-
{
92-
static_assert(TensorTools::TensorTraits<T>::rank == TensorTools::TensorTraits<T2>::rank,
93-
"Mathematical types must be same for arguments to relativelyFuzzEqual");
94-
if constexpr (TensorTools::TensorTraits<T>::rank == 0)
95-
return absolute_fuzzy_equals(
96-
var1,
97-
var2,
98-
tol * (std::abs(MetaPhysicL::raw_value(var1)) + std::abs(MetaPhysicL::raw_value(var2))));
99-
else if constexpr (TensorTools::TensorTraits<T>::rank == 1)
100-
{
101-
for (const auto i : make_range(libmesh_dim))
102-
if (!relative_fuzzy_equals(var1(i), var2(i), tol))
103-
return false;
104-
105-
return true;
106-
}
107-
else if constexpr (TensorTools::TensorTraits<T>::rank == 2)
108-
{
109-
for (const auto i : make_range(libmesh_dim))
110-
for (const auto j : make_range(libmesh_dim))
111-
if (!relative_fuzzy_equals(var1(i, j), var2(i, j), tol))
112-
return false;
113-
114-
return true;
115-
}
116-
}
117-
else
118-
{
119-
// We dare to dream
120-
libmesh_assert(var1.size() == var2.size());
121-
for (const auto i : index_range(var1))
122-
if (!relative_fuzzy_equals(var1(i), var2(i), tol))
123-
return false;
124-
125-
return true;
126-
}
94+
return absolute_fuzzy_equals(
95+
var1,
96+
var2,
97+
tol * (MetaPhysicL::raw_value(l1_norm(var1)) + MetaPhysicL::raw_value(l1_norm(var2))));
12798
}
128-
12999
#else
130-
131-
template <typename T,
132-
typename T2,
133-
typename std::enable_if<ScalarTraits<T>::value && ScalarTraits<T2>::value, int>::type>
100+
/**
101+
* Function to check whether two variables are equal within an absolute tolerance
102+
* @param var1 The first variable to be checked
103+
* @param var2 The second variable to be checked
104+
* @param tol The tolerance to be used
105+
* @return true if var1 and var2 are equal within tol
106+
*/
107+
template <typename T, typename T2>
134108
bool
135-
absolute_fuzzy_equals(const T & var1, const T2 & var2, const Real tol)
109+
absolute_fuzzy_equals(const T & var1, const T2 & var2, const Real tol = TOLERANCE * TOLERANCE)
136110
{
137-
return (std::abs(var1 - var2) <= tol);
111+
return l1_norm_diff(var1, var2) <= tol;
138112
}
139113

140114
/**
@@ -148,43 +122,9 @@ template <typename T, typename T2>
148122
bool
149123
relative_fuzzy_equals(const T & var1, const T2 & var2, const Real tol = TOLERANCE * TOLERANCE)
150124
{
151-
if constexpr (ScalarTraits<T>::value || TensorTools::MathWrapperTraits<T>::value)
152-
{
153-
static_assert(TensorTools::TensorTraits<T>::rank == TensorTools::TensorTraits<T2>::rank,
154-
"Mathematical types must be same for arguments to relativelyFuzzEqual");
155-
if constexpr (TensorTools::TensorTraits<T>::rank == 0)
156-
return absolute_fuzzy_equals(var1, var2, tol * (std::abs(var1) + std::abs(var2)));
157-
else if constexpr (TensorTools::TensorTraits<T>::rank == 1)
158-
{
159-
for (const auto i : make_range(libmesh_dim))
160-
if (!relative_fuzzy_equals(var1(i), var2(i), tol))
161-
return false;
162-
163-
return true;
164-
}
165-
else if constexpr (TensorTools::TensorTraits<T>::rank == 2)
166-
{
167-
for (const auto i : make_range(libmesh_dim))
168-
for (const auto j : make_range(libmesh_dim))
169-
if (!relative_fuzzy_equals(var1(i, j), var2(i, j), tol))
170-
return false;
171-
172-
return true;
173-
}
174-
}
175-
else
176-
{
177-
// We dare to dream
178-
mooseAssert(var1.size() == var2.size(), "These must be the same size");
179-
for (const auto i : index_range(var1))
180-
if (!relative_fuzzy_equals(var1(i), var2(i), tol))
181-
return false;
182-
183-
return true;
184-
}
125+
return absolute_fuzzy_equals(var1, var2, tol * (l1_norm(var1) + l1_norm(var2)));
185126
}
186-
187-
#endif // !LIBMESH_HAVE_METAPHYSICL
127+
#endif
188128

189129
} // namespace libMesh
190130

tests/numerics/type_vector_test.h

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -2,6 +2,7 @@
22
#define TYPE_VECTOR_TEST_H
33

44
#include <libmesh/type_vector.h>
5+
#include <libmesh/fuzzy_equals.h>
56

67
#include "libmesh_cppunit.h"
78

@@ -364,6 +365,7 @@ class TypeVectorTestBase : public CppUnit::TestCase {
364365

365366
CPPUNIT_ASSERT( (*basem_1_1_1) == (*basem_1_1_1) );
366367
CPPUNIT_ASSERT( !((*basem_1_1_1) == (*basem_n1_1_n1)) );
368+
CPPUNIT_ASSERT( relative_fuzzy_equals(*basem_1_1_1, *basem_1_1_1) );
367369
}
368370

369371
void testInEqualityBase()
@@ -372,6 +374,7 @@ class TypeVectorTestBase : public CppUnit::TestCase {
372374

373375
CPPUNIT_ASSERT( !((*basem_1_1_1) != (*basem_1_1_1)) );
374376
CPPUNIT_ASSERT( (*basem_1_1_1) != (*basem_n1_1_n1) );
377+
CPPUNIT_ASSERT( !relative_fuzzy_equals(*basem_1_1_1, *basem_n1_1_n1) );
375378
}
376379

377380
void testAssignmentBase()

0 commit comments

Comments
 (0)