From bd013067eca331fe53bfbf1ce8e32dc301ef172f Mon Sep 17 00:00:00 2001 From: Oliver Sander <oliver.sander@tu-dresden.de> Date: Tue, 14 Apr 2020 20:58:16 +0200 Subject: [PATCH] Fix adolcmaterialtest.cc The test erroneously checked for a small *absolute* error, but the numbers in this test are very large (around 1e8), so checking for a small *relative* error is more appropriate. This change makes the test pass. --- test/adolcmaterialtest.cc | 30 ++++++++++++++++++------------ 1 file changed, 18 insertions(+), 12 deletions(-) diff --git a/test/adolcmaterialtest.cc b/test/adolcmaterialtest.cc index 01344a0..3066260 100644 --- a/test/adolcmaterialtest.cc +++ b/test/adolcmaterialtest.cc @@ -1,6 +1,5 @@ #include<config.h> -#include <fstream> #include <random> #include <dune/fufem/assemblers/localassemblers/adolclocalenergy.hh> @@ -11,7 +10,6 @@ #include <dune/istl/bvector.hh> #include <dune/istl/bcrsmatrix.hh> -#include <dune/istl/matrixindexset.hh> #include <dune/grid/uggrid.hh> @@ -90,10 +88,15 @@ int main (int argc, char *argv[]) try for (size_t i=0; i < adolcGradient.size(); ++i) { - auto diff = adolcGradient[i] - paperGrad[i]; + auto relativeError = (adolcGradient[i] - paperGrad[i]).two_norm() + / std::max(adolcGradient[i].two_norm(), paperGrad[i].two_norm()); - if (diff.two_norm()>1e-5) - DUNE_THROW(Dune::Exception,"Wrong local derivative, error is "<<diff.two_norm()); + if (relativeError > 1e-11) + { + std::cout << "adolc:\n" << adolcGradient[i] << std::endl; + std::cout << "closed-form:\n" << paperGrad[i] << std::endl; + DUNE_THROW(Dune::Exception,"Wrong local derivative, error is " << relativeError); + } } // Test hessian assembler @@ -125,7 +128,6 @@ int main (int argc, char *argv[]) try auto adolcEndIt = adolcRow.end(); auto paperIt = paperRow.begin(); - auto paperEndIt = paperRow.end(); for (; adolcIt != adolcEndIt; ++adolcIt, ++paperIt) { @@ -133,18 +135,22 @@ int main (int argc, char *argv[]) try DUNE_THROW(Dune::Exception,"Not the same sparsity pattern!"<<adolcIt.index() <<"!="<<paperIt.index()); - auto diff = *adolcIt; - diff -= *paperIt; + auto relativeError = (*adolcIt - *paperIt).frobenius_norm() + / std::max(adolcIt->frobenius_norm(), paperIt->frobenius_norm()); - if (diff.frobenius_norm() > 1e-4) - DUNE_THROW(Dune::Exception,"Wrong local hessian, error is "<<diff.frobenius_norm()); + if (relativeError > 1e-11) + { + std::cout << "adolc:\n" << *adolcIt << std::endl; + std::cout << "closed-form:\n" << *paperIt << std::endl; + DUNE_THROW(Dune::Exception,"Wrong local hessian, error is " << relativeError); + } } - assert(paperIt==paperEndIt); + assert(paperIt==paperRow.end()); } return 0; -} catch (Exception e) { +} catch (const Exception& e) { std::cout << e << std::endl; return 1; } -- GitLab