Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Stepper_ErrorNorm_impl.hpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef Tempus_Stepper_ErrorNorm_impl_hpp
11 #define Tempus_Stepper_ErrorNorm_impl_hpp
12 
14 #include "Thyra_DefaultSerialDenseLinearOpWithSolve_decl.hpp"
15 #include "Thyra_MultiVectorStdOps_decl.hpp"
16 #include "Thyra_VectorSpaceBase_decl.hpp"
17 #include "Thyra_VectorStdOps_decl.hpp"
18 
20 
21 namespace Tempus {
22 
23 template <class Scalar>
25  : relTol_(1.0e-4), abssTol_(1.0e-4)
26 {
27 }
28 
29 template <class Scalar>
31  const Scalar absTol)
32  : relTol_(relTol), abssTol_(absTol)
33 {
34 }
35 
36 template <class Scalar>
39  const Teuchos::RCP<const Thyra::VectorBase<Scalar>> &xNext,
40  const Teuchos::RCP<const Thyra::VectorBase<Scalar>> &err)
41 {
42  if (errorWeightVector_ == Teuchos::null)
43  errorWeightVector_ = Thyra::createMember(x->space());
44 
45  if (u_ == Teuchos::null) u_ = Thyra::createMember(x->space());
46 
47  if (uNext_ == Teuchos::null) uNext_ = Thyra::createMember(x->space());
48 
49  // Compute: Atol + max(|u^n|, |u^{n+1}| ) * Rtol
50  Thyra::abs(*x, u_.ptr());
51  Thyra::abs(*xNext, uNext_.ptr());
52  Thyra::pair_wise_max_update(relTol_, *u_, uNext_.ptr());
53  if (!approxZero(abssTol_)) {
54  Thyra::add_scalar(abssTol_, uNext_.ptr());
55  }
56  else {
57  Scalar absTol = Thyra::norm_2(*uNext_) * numericalTol<Scalar>();
58  if (approxZero(absTol)) absTol = numericalTol<Scalar>();
59  Thyra::add_scalar(absTol, uNext_.ptr());
60  }
61 
62  Thyra::assign(errorWeightVector_.ptr(),
64  Thyra::ele_wise_divide(Teuchos::as<Scalar>(1.0), *err, *uNext_,
65  errorWeightVector_.ptr());
66 
67  const auto space_dim = err->space()->dim();
68  Scalar err_norm = std::abs(Thyra::norm(*errorWeightVector_) / space_dim);
69  return err_norm;
70 }
71 
72 template <class Scalar>
75 {
76  if (scratchVector_ == Teuchos::null)
77  scratchVector_ = Thyra::createMember(x->space());
78 
79  Thyra::assign(scratchVector_.ptr(), *x); // | U |
80  Thyra::abs(*x, scratchVector_.ptr());
81  Thyra::Vt_S(scratchVector_.ptr(), relTol_);
82  if (!approxZero(abssTol_)) {
83  Thyra::Vp_S(scratchVector_.ptr(), abssTol_);
84  }
85  else {
86  Scalar absTol = Thyra::norm_2(*scratchVector_) * numericalTol<Scalar>();
87  if (approxZero(absTol)) absTol = numericalTol<Scalar>();
88  Thyra::add_scalar(absTol, scratchVector_.ptr());
89  }
90  Thyra::ele_wise_divide(Teuchos::as<Scalar>(1.0), *x, *scratchVector_,
91  scratchVector_.ptr());
92  Scalar err = Thyra::norm_inf(*scratchVector_);
93  return err;
94 }
95 
96 } // namespace Tempus
97 
98 #endif
Scalar computeWRMSNorm(const Teuchos::RCP< const Thyra::VectorBase< Scalar >> &x, const Teuchos::RCP< const Thyra::VectorBase< Scalar >> &xNext, const Teuchos::RCP< const Thyra::VectorBase< Scalar >> &err)
Compute the weigthed root mean square norm.
bool approxZero(Scalar value, Scalar tol=Teuchos::ScalarTraits< Scalar >::sfmin())
Test if value is approximately zero within tolerance.
Scalar errorNorm(const Teuchos::RCP< const Thyra::VectorBase< Scalar >> &x)
Compute the error Norm.