Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
02_Use_ModelEvaluator.cpp
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#include <iomanip>
11#include <iostream>
12#include <stdlib.h>
13#include <math.h>
14#include "Teuchos_StandardCatchMacros.hpp"
15
16#include "Thyra_VectorStdOps.hpp"
17#include "Thyra_DefaultSpmdVectorSpace.hpp"
18#include "Thyra_DetachedVectorView.hpp"
19
21
22
23using namespace std;
24using Teuchos::RCP;
25
67int main(int argc, char *argv[])
68{
69 bool verbose = true;
70 bool success = false;
71 try {
72 // Construct ModelEvaluator
73 Teuchos::RCP<const Thyra::ModelEvaluator<double> >
74 model = Teuchos::rcp(new VanDerPol_ModelEvaluator_02<double>());
75
76 // Get initial conditions from ModelEvaluator::getNominalValues.
77 int n = 0;
78 double time = 0.0;
79 bool passed = true; // ICs are considered passed.
80 RCP<Thyra::VectorBase<double> > x_n =
81 model->getNominalValues().get_x()->clone_v();
82 RCP<Thyra::VectorBase<double> > xDot_n =
83 model->getNominalValues().get_x_dot()->clone_v();
84
85 // Timestep size
86 double finalTime = 2.0;
87 int nTimeSteps = 2001;
88 const double constDT = finalTime/(nTimeSteps-1);
89
90 // Advance the solution to the next timestep.
91 cout << n << " " << time << " " << get_ele(*(x_n), 0)
92 << " " << get_ele(*(x_n), 1) << endl;
93 while (passed && time < finalTime && n < nTimeSteps) {
94
95 // Initialize next time step
96 RCP<Thyra::VectorBase<double> > x_np1 = x_n->clone_v(); // at time index n+1
97
98 // Set the timestep and time.
99 double dt = constDT;
100 time = (n+1)*dt;
101
102 // For explicit ODE formulation, xDot = f(x, t),
103 // xDot is part of the outArgs.
104 auto inArgs = model->createInArgs();
105 auto outArgs = model->createOutArgs();
106 inArgs.set_t(time);
107 inArgs.set_x(x_n);
108 inArgs.set_x_dot(Teuchos::null);
109 outArgs.set_f(xDot_n);
110
111 // Righthand side evaluation and time-derivative at n.
112 model->evalModel(inArgs, outArgs);
113
114 // Take the timestep - Forward Euler
115 Thyra::V_VpStV(x_np1.ptr(), *x_n, dt, *xDot_n);
116
117 // Test if solution has passed.
118 if ( std::isnan(Thyra::norm(*x_np1)) ) {
119 passed = false;
120 } else {
121 // Promote to next step (n <- n+1).
122 Thyra::V_V(x_n.ptr(), *x_np1);
123 n++;
124 }
125
126 // Output
127 if ( n%100 == 0 )
128 cout << n << " " << time << " " << get_ele(*(x_n), 0)
129 << " " << get_ele(*(x_n), 1) << endl;
130 }
131
132 // Test for regression.
133 RCP<Thyra::VectorBase<double> > x_regress = x_n->clone_v();
134 {
135 Thyra::DetachedVectorView<double> x_regress_view(*x_regress);
136 x_regress_view[0] = -1.59496108218721311;
137 x_regress_view[1] = 0.96359412806611255;
138 }
139
140 RCP<Thyra::VectorBase<double> > x_error = x_n->clone_v();
141 Thyra::V_VmV(x_error.ptr(), *x_n, *x_regress);
142 double x_L2norm_error = Thyra::norm_2(*x_error );
143 double x_L2norm_regress = Thyra::norm_2(*x_regress);
144
145 cout << "Relative L2 Norm of the error (regression) = "
146 << x_L2norm_error/x_L2norm_regress << endl;
147 if ( x_L2norm_error > 1.0e-08*x_L2norm_regress) {
148 passed = false;
149 cout << "FAILED regression constraint!" << endl;
150 }
151
152 RCP<Thyra::VectorBase<double> > x_best = x_n->clone_v();
153 {
154 Thyra::DetachedVectorView<double> x_best_view(*x_best);
155 x_best_view[0] = -1.59496108218721311;
156 x_best_view[1] = 0.96359412806611255;
157 }
158
159 Thyra::V_VmV(x_error.ptr(), *x_n, *x_best);
160 x_L2norm_error = Thyra::norm_2(*x_error);
161 double x_L2norm_best = Thyra::norm_2(*x_best );
162
163 cout << "Relative L2 Norm of the error (best) = "
164 << x_L2norm_error/x_L2norm_best << endl;
165 if ( x_L2norm_error > 0.02*x_L2norm_best) {
166 passed = false;
167 cout << "FAILED best constraint!" << endl;
168 }
169 if (passed) success = true;
170 }
171 TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
172
173 if(success)
174 cout << "\nEnd Result: Test Passed!" << std::endl;
175
176 return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
177}
int main(int argc, char *argv[])
ModelEvaluator implementation for the example van der Pol Problem.