Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
taylor_example.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Sacado Package
4 //
5 // Copyright 2006 NTESS and the Sacado contributors.
6 // SPDX-License-Identifier: LGPL-2.1-or-later
7 // *****************************************************************************
8 // @HEADER
9 
10 // taylor_example
11 //
12 // usage:
13 // taylor_example
14 //
15 // output:
16 // prints the results of computing a single Taylor series expansion of
17 // the solution to:
18 //
19 // dx/dt = 1 + x^2, x(0) = 1.0;
20 //
21 // The exact solution is x(t) = tan(t + pi/4)
22 
23 #include <iostream>
24 
25 #include "Sacado_No_Kokkos.hpp"
26 
27 // Function implementing RHS of ODE
28 template <typename ScalarT>
29 void func(ScalarT& f, const ScalarT& x) {
30  f = 1.0 + x*x;
31 }
32 
33 int main(int argc, char **argv)
34 {
35  double x0 = 1.0; // Initial condition
36  int deg = 40; // Degree of Taylor series solution
37 
38  Sacado::Tay::Taylor<double> x = x0; // Taylor polynomial for independent
39  Sacado::Tay::Taylor<double> f; // Taylor polynomial for dependent
40 
41  // Reserve space for degree deg coefficients
42  x.reserve(deg);
43 
44  // Compute Taylor series solution to dx/dt = f(x)
45  for (int k=0; k<deg; k++) {
46  func(f, x);
47 
48  // Set next coefficient
49  x.resize(k+1, true);
50 
51  // x_{k+1} = f_k / (k+1)
52  x.fastAccessCoeff(k+1) = f.coeff(k) / (k+1);
53  }
54 
55  // Print Taylor series solution
56  std::cout << "Taylor series solution = " << std::endl
57  << x << std::endl;
58 
59  // Compute Taylor series expansion of solution x(t) = tan(t+pi/4)
60  double pi = std::atan(1.0)*4.0;
62  t.fastAccessCoeff(0) = pi/4.0;
63  t.fastAccessCoeff(1) = 1.0;
65 
66  // Print expansion of solution
67  std::cout << "Exact expansion = " << std::endl
68  << u << std::endl;
69 
70  // Compute maximum relative error
71  double max_err = 0.0;
72  double err = 0.0;
73  for (int k=0; k<=deg; k++) {
74  err = std::fabs(x.coeff(k) - u.coeff(k)) / (1.0 + fabs(u.coeff(k)));
75  if (err > max_err) max_err = err;
76  }
77  std::cout << "Maximum relative error = " << max_err << std::endl;
78 
79  double tol = 1.0e-12;
80  if (max_err < tol){
81  std::cout << "\nExample passed!" << std::endl;
82  return 0;
83  }
84  else {
85  std::cout <<"\nSomething is wrong, example failed!" << std::endl;
86  return 1;
87  }
88 }
89 
90 
91 
void resize(int d, bool keep_coeffs)
Resize polynomial to degree d.
void f()
atan(expr.val())
T & fastAccessCoeff(int i)
Returns degree i term without bounds checking.
int main()
Definition: ad_example.cpp:171
tan(expr.val())
const T * coeff() const
Returns Taylor coefficient array.
void reserve(int d)
Reserve space for a degree d polynomial.
const double tol
const T func(int n, T *x)
Definition: ad_example.cpp:29
fabs(expr.val())