Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
03_Intro_SolutionState.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
20#include "../02_Use_ModelEvaluator/VanDerPol_ModelEvaluator_02.hpp"
21
22#include "Tempus_SolutionState.hpp"
23
24
25using namespace std;
26using Teuchos::RCP;
27
75int main(int argc, char *argv[])
76{
77 bool verbose = true;
78 bool success = false;
79 try {
80 // Construct ModelEvaluator
81 Teuchos::RCP<const Thyra::ModelEvaluator<double> >
82 model = Teuchos::rcp(new VanDerPol_ModelEvaluator_02<double>());
83
84 // Setup initial condition SolutionState --------------------
85 auto solState = Tempus::createSolutionStateX(
86 model->getNominalValues().get_x()->clone_v());
87 solState->setIndex (0);
88 solState->setTime (0.0);
89 solState->setTimeStep(0.0); // By convention, the IC has dt=0.
90 solState->setSolutionStatus(Tempus::Status::PASSED); // ICs are considered passed.
91 RCP<Thyra::VectorBase<double> > xDot_n =
92 model->getNominalValues().get_x_dot()->clone_v();
93
94
95 // Timestep size
96 double finalTime = 2.0;
97 int nTimeSteps = 2001;
98 const double constDT = finalTime/(nTimeSteps-1);
99
100 // Advance the solution to the next timestep.
101 while (solState->getSolutionStatus() == Tempus::Status::PASSED &&
102 solState->getTime() < finalTime &&
103 solState->getIndex() < nTimeSteps) {
104
105 // Initialize next time step
106 RCP<Thyra::VectorBase<double> > x_n = solState->getX();
107 RCP<Thyra::VectorBase<double> > x_np1 = solState->getX()->clone_v(); // at time index n+1
108 solState->setSolutionStatus(Tempus::Status::WORKING);
109
110 // Set the timestep and time for the working solution i.e., n+1.
111 int index = solState->getIndex()+1;
112 double dt = constDT;
113 double time = index*dt;
114
115 // For explicit ODE formulation, xDot = f(x, t),
116 // xDot is part of the outArgs.
117 auto inArgs = model->createInArgs();
118 auto outArgs = model->createOutArgs();
119 inArgs.set_t(time);
120 inArgs.set_x(x_n);
121 inArgs.set_x_dot(Teuchos::null);
122 outArgs.set_f(xDot_n);
123
124 // Righthand side evaluation and time-derivative at n.
125 model->evalModel(inArgs, outArgs);
126
127 // Take the timestep - Forward Euler
128 Thyra::V_VpStV(x_np1.ptr(), *x_n, dt, *xDot_n);
129
130 // Test if solution has passed.
131 if ( std::isnan(Thyra::norm(*x_np1)) ) {
132 solState->setSolutionStatus(Tempus::Status::FAILED);
133 } else {
134 // Promote to next step (n <- n+1).
135 Thyra::V_V(x_n.ptr(), *x_np1);
136 solState->setIndex (index);
137 solState->setTime (time);
138 solState->setTimeStep(constDT);
139 solState->setSolutionStatus(Tempus::Status::PASSED);
140 }
141
142 // Output
143 if ( solState->getIndex()%100 == 0 )
144 cout << solState->getIndex() << " " << time
145 << " " << get_ele(*(x_n), 0)
146 << " " << get_ele(*(x_n), 1) << endl;
147 }
148
149 // Test for regression.
150 RCP<Thyra::VectorBase<double> > x_n = solState->getX();
151 RCP<Thyra::VectorBase<double> > x_regress = x_n->clone_v();
152 {
153 Thyra::DetachedVectorView<double> x_regress_view(*x_regress);
154 x_regress_view[0] = -1.59496108218721311;
155 x_regress_view[1] = 0.96359412806611255;
156 }
157
158 RCP<Thyra::VectorBase<double> > x_error = x_n->clone_v();
159 Thyra::V_VmV(x_error.ptr(), *x_n, *x_regress);
160 double x_L2norm_error = Thyra::norm_2(*x_error );
161 double x_L2norm_regress = Thyra::norm_2(*x_regress);
162
163 cout << "Relative L2 Norm of the error (regression) = "
164 << x_L2norm_error/x_L2norm_regress << endl;
165 if ( x_L2norm_error > 1.0e-08*x_L2norm_regress) {
166 solState->setSolutionStatus(Tempus::Status::FAILED);
167 cout << "FAILED regression constraint!" << endl;
168 }
169
170 RCP<Thyra::VectorBase<double> > x_best = x_n->clone_v();
171 {
172 Thyra::DetachedVectorView<double> x_best_view(*x_best);
173 x_best_view[0] = -1.59496108218721311;
174 x_best_view[1] = 0.96359412806611255;
175 }
176
177 Thyra::V_VmV(x_error.ptr(), *x_n, *x_best);
178 x_L2norm_error = Thyra::norm_2(*x_error);
179 double x_L2norm_best = Thyra::norm_2(*x_best );
180
181 cout << "Relative L2 Norm of the error (best) = "
182 << x_L2norm_error/x_L2norm_best << endl;
183 if ( x_L2norm_error > 0.02*x_L2norm_best) {
184 solState->setSolutionStatus(Tempus::Status::FAILED);
185 cout << "FAILED best constraint!" << endl;
186 }
187 if (solState->getSolutionStatus() == Tempus::Status::PASSED) success = true;
188 }
189 TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
190
191 if(success)
192 cout << "\nEnd Result: Test Passed!" << std::endl;
193
194 return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
195}
int main(int argc, char *argv[])
ModelEvaluator implementation for the example van der Pol Problem.
Teuchos::RCP< SolutionState< Scalar > > createSolutionStateX(const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdot=Teuchos::null, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdotdot=Teuchos::null)
Nonmember constructor from non-const solution vectors, x.