$extrastylesheet
libMesh::TransientSystem< Base > Class Template Reference

#include <transient_system.h>

List of all members.

Public Types

typedef TransientSystem< Base > sys_type

Public Member Functions

 TransientSystem (EquationSystems &es, const std::string &name, const unsigned int number)
virtual ~TransientSystem ()
sys_typesystem ()
virtual void clear ()
virtual void reinit ()
virtual std::string system_type () const
Number old_solution (const dof_id_type global_dof_number) const
Number older_solution (const dof_id_type global_dof_number) const

Public Attributes

UniquePtr< NumericVector
< Number > > 
old_local_solution
UniquePtr< NumericVector
< Number > > 
older_local_solution

Protected Member Functions

virtual void init_data ()
virtual void re_update ()

Detailed Description

template<class Base>
class libMesh::TransientSystem< Base >

This class provides a specific system class. It aims at transient systems, offering nothing more than just the essentials needed to solve a system. Note that still additional vectors/matrices may be added, as offered in the parent classes.

Definition at line 48 of file transient_system.h.


Member Typedef Documentation

template<class Base>
typedef TransientSystem<Base> libMesh::TransientSystem< Base >::sys_type

The type of system.

Reimplemented in libMesh::TransientRBConstruction.

Definition at line 68 of file transient_system.h.


Constructor & Destructor Documentation

template<class Base >
libMesh::TransientSystem< Base >::TransientSystem ( EquationSystems es,
const std::string &  name,
const unsigned int  number 
)

Constructor. Initializes required data structures.

Definition at line 38 of file transient_system.C.

References libMesh::GHOSTED, libMesh::TransientSystem< Base >::old_local_solution, libMesh::TransientSystem< Base >::older_local_solution, and libMesh::SERIAL.

                                                                      :

  Base                 (es, name_in, number_in)
{
#ifdef LIBMESH_ENABLE_GHOSTED
  old_local_solution =
    UniquePtr<NumericVector<Number> >
    (&(this->add_vector("_transient_old_local_solution", true, GHOSTED)));
  older_local_solution =
    UniquePtr<NumericVector<Number> >
    (&(this->add_vector("_transient_older_local_solution", true, GHOSTED)));
#else
  old_local_solution =
    UniquePtr<NumericVector<Number> >
    (&(this->add_vector("_transient_old_local_solution", true, SERIAL)));
  older_local_solution =
    UniquePtr<NumericVector<Number> >
    (&(this->add_vector("_transient_older_local_solution", true, SERIAL)));
#endif
}
template<class Base >
libMesh::TransientSystem< Base >::~TransientSystem ( ) [virtual]

Destructor.

Definition at line 64 of file transient_system.C.

{
  this->clear();

  // We still have UniquePtrs for API compatibility, but
  // now that we're System::add_vector()ing these, we can trust
  // the base class to handle memory management
  old_local_solution.release();
  older_local_solution.release();
}

Member Function Documentation

template<class Base >
void libMesh::TransientSystem< Base >::clear ( ) [virtual]

Clear all the data structures associated with the system.

Reimplemented in libMesh::TransientRBConstruction.

Definition at line 78 of file transient_system.C.

{
  // clear the parent data
  Base::clear();

  // the old & older local solutions
  // are now deleted by System!
  // old_local_solution->clear();
  // older_local_solution->clear();

  // FIXME: This preserves maximum backwards compatibility,
  // but is probably grossly unnecessary:
  old_local_solution.release();
  older_local_solution.release();

  old_local_solution =
    UniquePtr<NumericVector<Number> >
    (&(this->add_vector("_transient_old_local_solution")));
  older_local_solution =
    UniquePtr<NumericVector<Number> >
    (&(this->add_vector("_transient_older_local_solution")));
}
template<class Base >
void libMesh::TransientSystem< Base >::init_data ( ) [protected, virtual]

Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.

Definition at line 105 of file transient_system.C.

References libMesh::get_dof_map(), libMesh::GHOSTED, and libMesh::SERIAL.

{
  // initialize parent data
  Base::init_data();

  // Initialize the old & older solutions
  // Using new ghosted vectors if enabled
#ifdef LIBMESH_ENABLE_GHOSTED
  old_local_solution->init   (this->n_dofs(), this->n_local_dofs(),
                              this->get_dof_map().get_send_list(), false,
                              GHOSTED);
  older_local_solution->init (this->n_dofs(), this->n_local_dofs(),
                              this->get_dof_map().get_send_list(), false,
                              GHOSTED);
#else
  old_local_solution->init   (this->n_dofs(), false, SERIAL);
  older_local_solution->init (this->n_dofs(), false, SERIAL);
#endif
}
template<class Base >
Number libMesh::TransientSystem< Base >::old_solution ( const dof_id_type  global_dof_number) const
Returns:
the old solution (at the previous timestep) for the specified global DOF.

Definition at line 177 of file transient_system.C.

References libMesh::get_dof_map().

{
  // Check the sizes
  libmesh_assert_less (global_dof_number, this->get_dof_map().n_dofs());
  libmesh_assert_less (global_dof_number, old_local_solution->size());

  return (*old_local_solution)(global_dof_number);
}
template<class Base >
Number libMesh::TransientSystem< Base >::older_solution ( const dof_id_type  global_dof_number) const
Returns:
the older solution (two timesteps ago) for the specified global DOF.

Definition at line 189 of file transient_system.C.

References libMesh::get_dof_map().

{
  // Check the sizes
  libmesh_assert_less (global_dof_number, this->get_dof_map().n_dofs());
  libmesh_assert_less (global_dof_number, older_local_solution->size());

  return (*older_local_solution)(global_dof_number);
}
template<class Base >
void libMesh::TransientSystem< Base >::re_update ( ) [protected, virtual]

Re-update the local values when the mesh has changed. This method takes the data updated by update() and makes it up-to-date on the current mesh.

Definition at line 142 of file transient_system.C.

References libMesh::get_dof_map().

{
  // re_update the parent system
  Base::re_update ();

  const std::vector<dof_id_type>& send_list = this->get_dof_map().get_send_list ();

  const dof_id_type first_local_dof = Base::get_dof_map().first_dof();
  const dof_id_type end_local_dof  = Base::get_dof_map().end_dof();

  // Check sizes
  libmesh_assert_greater_equal (end_local_dof, first_local_dof);
  libmesh_assert_greater_equal (older_local_solution->size(), send_list.size());
  libmesh_assert_greater_equal (old_local_solution->size(), send_list.size());

  // Even if we don't have to do anything ourselves, localize() may
  // use parallel_only tools
  // if (first_local_dof == end_local_dof)
  //   return;

  // Update the old & older solutions with the send_list,
  // which may have changed since their last update.
  older_local_solution->localize (first_local_dof,
                                  end_local_dof-1,
                                  send_list);

  old_local_solution->localize (first_local_dof,
                                end_local_dof-1,
                                send_list);
}
template<class Base >
void libMesh::TransientSystem< Base >::reinit ( ) [virtual]

Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be used.

Definition at line 128 of file transient_system.C.

{
  // initialize parent data
  Base::reinit();

  // Project the old & older vectors to the new mesh
  // The System::reinit handles this now
  // this->project_vector (*old_local_solution);
  // this->project_vector (*older_local_solution);
}
template<class Base>
sys_type& libMesh::TransientSystem< Base >::system ( ) [inline]
Returns:
a clever pointer to the system.

Definition at line 73 of file transient_system.h.

{ return *this; }
template<class Base >
std::string libMesh::TransientSystem< Base >::system_type ( ) const [inline, virtual]
Returns:
"Transient" prepended to T::system_type(). Helps in identifying the system type in an equation system file.

Definition at line 160 of file transient_system.h.

{
  std::string type = "Transient";
  type += Base::system_type ();

  return type;
}

Member Data Documentation

template<class Base>
UniquePtr<NumericVector<Number> > libMesh::TransientSystem< Base >::old_local_solution

All the values I need to compute my contribution to the simulation at hand. Think of this as the current solution with any ghost values needed from other processors.

Definition at line 116 of file transient_system.h.

Referenced by libMesh::TransientSystem< Base >::TransientSystem().

template<class Base>
UniquePtr<NumericVector<Number> > libMesh::TransientSystem< Base >::older_local_solution

All the values I need to compute my contribution to the simulation at hand. Think of this as the current solution with any ghost values needed from other processors.

Definition at line 124 of file transient_system.h.

Referenced by libMesh::TransientSystem< Base >::TransientSystem().


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