Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Transition from Example 0 to Example 1

Goal

The primary goal of Example 1: Utilize Thyra is to replace the raw arrays with Thyra vectors, while preserving the same van der Pol problem and the same Forward Euler stepping logic used in Example 0: Basic Problem.

New concepts introduced from Teuchos and Thyra are

  • Teuchos::RCP
    • Teuchos::RCP memory management can be compared to std::shared_ptr
  • Thyra::VectorBase
    • Replaces raw C++ arrays
  • Thyra::VectorSpaceBase
    • Introduced to define the state space.
  • Thyra::DetactedVectorView and Thyra::ConstDetactedVectorView
    • Detached vector views are used for direct element access.
  • Thyra::V_VpSV, Thra::_V, Thyra::norm, Thyra::norm_2
    • Thyra helper routines to express vector algebra.

This allows Tempus and other Trilinos packages to operate on abstract numerical interfaces rather than concrete array types. Moving from raw arrays to Thyra vectors is the first step toward interoperability with Tempus steppers, model evaluators, and solver infrastructure.

What stays the same

  • The van der Pol equations are unchanged.
  • The timestepper is still Forward Euler.
  • The time loop structure is the same.
  • The regression check follows the same pattern.

Selected code excerpts

Below are select code snippets and comments on the uses Teuchos::RCP and Thyra functionality.


Memory Allocation

Setup Thyra vectors. The raw double arrays are replaced by Thyra vectors, which are created from a Thyra::VectorSpaceBase. This step also introduces Teuchos::RCP, Trilinos's reference-counted smart pointer.

Before

double x_n[2];
double xDot_n[2];

After

int vectorLength = 2;
RCP<const Thyra::VectorSpaceBase<double> > xSpace =
Thyra::defaultSpmdVectorSpace<double>(vectorLength);
RCP<Thyra::VectorBase<double> > x_n = Thyra::createMember(xSpace);
RCP<Thyra::VectorBase<double> > xDot_n = Thyra::createMember(xSpace);

The state is now represented through abstract vector interfaces rather than fixed-size arrays.


Assign Initial Conditions

Initialize Thyra vectors. Initial values are assigned through Thyra::DetachedVectorView. The local scope ensures the detached views are destroyed immediately after use.

Before

x_n [0] = 2.0;
x_n [1] = 0.0;
xDot_n[0] = 0.0;
xDot_n[1] = -2.0/epsilon;

After

{ // Scope to delete DetachedVectorViews
Thyra::DetachedVectorView<double> x_n_view(*x_n);
x_n_view[0] = 2.0;
x_n_view[1] = 0.0;
Thyra::DetachedVectorView<double> xDot_n_view(*xDot_n);
xDot_n_view[0] = 0.0;
xDot_n_view[1] = -2.0/epsilon;
}

Evaluate the RHS

Evaluating the Right-Hand Side The van der Pol right-hand side is still evaluated directly in the application code, but element access now goes through detached vector views.

Before

xDot_n[0] = x_n[1];
xDot_n[1] = ((1.0 - x_n[0]*x_n[0])*x_n[1] - x_n[0])/epsilon;

After

{
Thyra::ConstDetachedVectorView<double> x_n_view(*x_n);
Thyra::DetachedVectorView<double> xDot_n_view(*xDot_n);
xDot_n_view[0] = x_n_view[1];
xDot_n_view[1] =
((1.0-x_n_view[0]*x_n_view[0])*x_n_view[1]-x_n_view[0])/epsilon;
}

Element-wise access is still possible, but it now uses Thyra utilities.


Perform Forward Euler

Use Thyra vector algebra. The Forward Euler update is expressed with Thyra::V_VpStV, which performs the vector operation $x^{n+1} = x^n + dt\,\dot{x}^n$.

Before

x_np1[0] = x_n[0] + dt*xDot_n[0];
x_np1[1] = x_n[1] + dt*xDot_n[1];

After

Thyra::V_VpStV(x_np1.ptr(), *x_n, dt, *xDot_n);

Promote to the next time step

Use Thyra vector assignment. The accepted solution is copied into the current state using Thyra::V_V.

Before

x_n[0] = x_np1[0];
x_n[1] = x_np1[1];
n++;

After

Thyra::V_V(x_n.ptr(), *x_np1);
n++;

How to compare the examples

A more detailed comparison can be viewed by diffing:

From the packages/tempus directory, a focused comparison of the main time-integration logic between these two examples can be generated locally in bash or zsh with:

git diff --no-index \
<(sed -n '81,141p' examples/00_Basic_Problem/00_Basic_Problem.cpp) \
<(sed -n '60,130p' examples/01_Utilize_Thyra/01_Utilize_Thyra.cpp)

This ignores leading header lines (for example, #include statements and Doxygen comments) and trailing regression-testing lines.

← Previous Example | Current Example | Next Example →