Doxygen API reference documentation for ideal.II
Public Member Functions | List of all members
idealii::TimeIteratorCollection< dim > Class Template Reference

A collection of slab iterators to simplify time marching. More...

#include <time_iterator.hh>

Public Member Functions

 TimeIteratorCollection ()
 Default constructor. More...
 
void add_iterator (slab::TriaIterator< dim > *it, spacetime::Triangulation< dim > *collection)
 Add a slab::TriaIterator<dim> iterator. More...
 
void add_iterator (slab::parallel::distributed::TriaIterator< dim > *it, spacetime::parallel::distributed::Triangulation< dim > *collection)
 Add a parallel::distributed::slab::TriaIterator<dim> iterator. More...
 
void add_iterator (slab::DoFHandlerIterator< dim > *it, spacetime::DoFHandler< dim > *collection)
 Add a slab::DoFHandlerIterator<dim> iterator. More...
 
void add_iterator (slab::VectorIterator< double > *it, spacetime::Vector< double > *collection)
 Add a slab::VectorIterator<dim> iterator. More...
 
void add_iterator (slab::TrilinosVectorIterator *it, spacetime::TrilinosVector *collection)
 Add a slab::TrilinosVectorIterator<dim> iterator. More...
 
void increment ()
 
void decrement ()
 
bool at_end ()
 
bool before_begin ()
 

Detailed Description

template<int dim>
class idealii::TimeIteratorCollection< dim >

A collection of slab iterators to simplify time marching.

The second type of spacetime objects like spacetime::Triangulation, spacetime::DoFHandler and spacetime::Vector are simply doubly linked lists of objects corresponding to a specific slab.

For forward or backward time marching these lists need to be traversed simultaneously. This class combines all iterator increments or decrements into a single function call.

Constructor & Destructor Documentation

◆ TimeIteratorCollection()

Default constructor.

Note
Any iterators have to be registered using one of the add_iterator functions.

Member Function Documentation

◆ add_iterator() [1/5]

template<int dim>
void idealii::TimeIteratorCollection< dim >::add_iterator ( slab::DoFHandlerIterator< dim > *  it,
spacetime::DoFHandler< dim > *  collection 
)

Add a slab::DoFHandlerIterator<dim> iterator.

Parameters
itA pointer to the slab::DoFHandlerIterator<dim>
collectionA pointer to the spacetime::DoFHandler<dim> whose member list includes it

◆ add_iterator() [2/5]

template<int dim>
void idealii::TimeIteratorCollection< dim >::add_iterator ( slab::parallel::distributed::TriaIterator< dim > *  it,
spacetime::parallel::distributed::Triangulation< dim > *  collection 
)

Add a parallel::distributed::slab::TriaIterator<dim> iterator.

Parameters
itA pointer to the parallel::distributed::slab::TriaIterator<dim>
collectionA pointer to the parallel::distributed::spacetime::Triangulation<dim> whose member list includes it

◆ add_iterator() [3/5]

template<int dim>
void idealii::TimeIteratorCollection< dim >::add_iterator ( slab::TriaIterator< dim > *  it,
spacetime::Triangulation< dim > *  collection 
)

Add a slab::TriaIterator<dim> iterator.

Parameters
itA pointer to the slab::TriaIterator<dim>
collectionA pointer to the spacetime::Triangulation<dim> whose member list includes it

◆ add_iterator() [4/5]

template<int dim>
void idealii::TimeIteratorCollection< dim >::add_iterator ( slab::TrilinosVectorIterator it,
spacetime::TrilinosVector collection 
)

Add a slab::TrilinosVectorIterator<dim> iterator.

Parameters
itA pointer to the slab::TrilinosVectorIterator<dim>
collectionA pointer to the spacetime::TrilinosVector<dim> whose member list includes it

◆ add_iterator() [5/5]

template<int dim>
void idealii::TimeIteratorCollection< dim >::add_iterator ( slab::VectorIterator< double > *  it,
spacetime::Vector< double > *  collection 
)

Add a slab::VectorIterator<dim> iterator.

Parameters
itA pointer to the slab::VectorIterator<dim>
collectionA pointer to the spacetime::Vector<dim> whose member list includes it

◆ at_end()

template<int dim>
bool idealii::TimeIteratorCollection< dim >::at_end ( )

For use as a stopping criterion in forward time marching.

Returns
true: if at least one of the iterators is a the corresponding end()
false: if none of the iterators is at the corresponding end()

◆ before_begin()

template<int dim>
bool idealii::TimeIteratorCollection< dim >::before_begin ( )

For use as a stopping criterion in backward time marching.

Returns
true: if at least one of the iterators is before the corresponding begin()
false: if none of the iterators is before the corresponding begin()

◆ decrement()

template<int dim>
void idealii::TimeIteratorCollection< dim >::decrement ( )

Decrements all added iterators via their prefix decrement operators.

◆ increment()

template<int dim>
void idealii::TimeIteratorCollection< dim >::increment ( )

Increments all added iterators via their prefix increment operators.


The documentation for this class was generated from the following file: