EquationSystems Class Reference

#include <equation_systems.h>

List of all members.

Public Types

enum  ReadFlags {
  READ_HEADER = 1, READ_DATA = 2, READ_ADDITIONAL_DATA = 4, READ_LEGACY_FORMAT = 8,
  TRY_READ_IFEMS = 16
}
enum  WriteFlags { WRITE_DATA = 1, WRITE_ADDITIONAL_DATA = 2, WRITE_PARALLEL_FILES = 4 }

Public Member Functions

 EquationSystems (MeshBase &mesh, MeshData *mesh_data=NULL)
virtual ~EquationSystems ()
void clear ()
void init ()
void reinit ()
void update ()
unsigned int n_systems () const
bool has_system (const std::string &name) const
template<typename T_sys >
const T_sys & get_system (const std::string &name) const
template<typename T_sys >
T_sys & get_system (const std::string &name)
template<typename T_sys >
const T_sys & get_system (const unsigned int num) const
template<typename T_sys >
T_sys & get_system (const unsigned int num)
const Systemget_system (const std::string &name) const
Systemget_system (const std::string &name)
const Systemget_system (const unsigned int num) const
Systemget_system (const unsigned int num)
Systemadd_system (const std::string &system_type, const std::string &name)
template<typename T_sys >
T_sys & add_system (const std::string &name)
void delete_system (const std::string &name)
unsigned int n_vars () const
unsigned int n_dofs () const
unsigned int n_active_dofs () const
virtual void solve ()
virtual void adjoint_solve (const QoISet &qoi_indices=QoISet())
virtual void sensitivity_solve (const ParameterVector &parameters)
void build_variable_names (std::vector< std::string > &var_names) const
void build_solution_vector (std::vector< Number > &soln, const std::string &system_name, const std::string &variable_name="all_vars") const
void build_solution_vector (std::vector< Number > &soln) const
void build_discontinuous_solution_vector (std::vector< Number > &soln) const
void read (const std::string &name, const libMeshEnums::XdrMODE, const unsigned int read_flags=(READ_HEADER|READ_DATA))
void write (const std::string &name, const libMeshEnums::XdrMODE, const unsigned int write_flags=(WRITE_DATA)) const
bool compare (const EquationSystems &other_es, const Real threshold, const bool verbose) const
std::string get_info () const
void print_info (std::ostream &os=std::cout) const
const MeshBaseget_mesh () const
MeshBaseget_mesh ()
bool has_mesh_data () const
const MeshDataget_mesh_data () const
MeshDataget_mesh_data ()
void allgather ()

Public Attributes

Parameters parameters

Protected Types

typedef std::map< std::string,
System * >::iterator 
system_iterator
typedef std::map< std::string,
System * >::const_iterator 
const_system_iterator

Protected Attributes

MeshBase_mesh
MeshData_mesh_data
std::map< std::string, System * > _systems

Private Member Functions

void _read_impl (const std::string &name, const libMeshEnums::XdrMODE, const unsigned int read_flags)
void _add_system_to_nodes_and_elems ()

Friends

std::ostream & operator<< (std::ostream &os, const EquationSystems &es)


Detailed Description

This is the EquationSystems class. It is in charge of handling all the various equation systems defined for a MeshBase. It may have multiple systems, which may be active or inactive, so that at different solution stages only a sub-set may be solved for. Also, through the templated access, different types of systems may be handled. Also other features, like flags, parameters, I/O etc are provided.

Author:
Benjamin S. Kirk, 2002-2007

Definition at line 65 of file equation_systems.h.


Member Typedef Documentation

typedef std::map<std::string, System*>::const_iterator EquationSystems::const_system_iterator [protected]

Typedef for constatnt system iterators

Definition at line 406 of file equation_systems.h.

typedef std::map<std::string, System*>::iterator EquationSystems::system_iterator [protected]

Typedef for system iterators

Definition at line 401 of file equation_systems.h.


Member Enumeration Documentation

Define enumeration to set properties in EquationSystems::read()

Enumerator:
READ_HEADER 
READ_DATA 
READ_ADDITIONAL_DATA 
READ_LEGACY_FORMAT 
TRY_READ_IFEMS 

Definition at line 72 of file equation_systems.h.

00072                  { READ_HEADER           = 1,
00073                    READ_DATA             = 2,
00074                    READ_ADDITIONAL_DATA  = 4,
00075                    READ_LEGACY_FORMAT    = 8,
00076                    TRY_READ_IFEMS        = 16 };

Define enumeration to set properties in EquationSystems::write()

Enumerator:
WRITE_DATA 
WRITE_ADDITIONAL_DATA 
WRITE_PARALLEL_FILES 

Definition at line 81 of file equation_systems.h.

00081                   { WRITE_DATA             = 1,
00082                     WRITE_ADDITIONAL_DATA  = 2,
00083                     WRITE_PARALLEL_FILES   = 4 };


Constructor & Destructor Documentation

EquationSystems::EquationSystems ( MeshBase mesh,
MeshData mesh_data = NULL 
)

Constructor.

Definition at line 50 of file equation_systems.C.

References parameters, and Parameters::set().

00050                                                                   :
00051   _mesh      (m),
00052   _mesh_data (mesh_data)
00053 {
00054   // Set default parameters
00055   this->parameters.set<Real>        ("linear solver tolerance") = TOLERANCE * TOLERANCE;
00056   this->parameters.set<unsigned int>("linear solver maximum iterations") = 5000;  
00057 }

EquationSystems::~EquationSystems (  )  [virtual]

Destructor. Should be virtual, since the user may want to derive subclasses of EquationSystems.

Definition at line 61 of file equation_systems.C.

References clear().

00062 {
00063   this->clear ();
00064 }


Member Function Documentation

void EquationSystems::_add_system_to_nodes_and_elems (  )  [private]

This function is used in the implementation of add_system, it loops over the nodes and elements of the Mesh, adding the system to each one. The main reason to separate this part is to avoid coupling this header file to mesh.h, and elem.h.

Definition at line 913 of file equation_systems.C.

References _mesh, MeshBase::elements_begin(), MeshBase::elements_end(), MeshBase::nodes_begin(), and MeshBase::nodes_end().

Referenced by add_system().

00914 {
00915   // All the nodes
00916   MeshBase::node_iterator       node_it  = _mesh.nodes_begin();
00917   const MeshBase::node_iterator node_end = _mesh.nodes_end();
00918  
00919   for ( ; node_it != node_end; ++node_it)
00920     (*node_it)->add_system();
00921  
00922   // All the elements
00923   MeshBase::element_iterator       elem_it  = _mesh.elements_begin();
00924   const MeshBase::element_iterator elem_end = _mesh.elements_end();
00925  
00926   for ( ; elem_it != elem_end; ++elem_it)
00927     (*elem_it)->add_system();
00928 }

void EquationSystems::_read_impl ( const std::string &  name,
const libMeshEnums::XdrMODE  mode,
const unsigned int  read_flags 
) [private]

Actual read implementation. This can be called repeatedly inside a try-catch block in an attempt to read broken files.

This program implements the output of an EquationSystems object. This warrants some documentation. The output file essentially consists of 11 sections:

     1.) A version header (for non-'legacy' formats, libMesh-0.7.0 and greater). 
     2.) The number of individual equation systems (unsigned int)
     
       for each system
                                                          
        3.)  The name of the system (string)
        4.)  The type of the system (string)
    
        handled through System::read():
    
     +-------------------------------------------------------------+
     |  5.) The number of variables in the system (unsigned int)   |
     |                                                             |
     |   for each variable in the system                           |
     |                                                             |
     |    6.) The name of the variable (string)                    |
     |                                                             |
     |    7.) Combined in an FEType:                               |
     |         - The approximation order(s) of the variable (Order |
     |           Enum, cast to int/s)                              |
     |         - The finite element family/ies of the variable     |
     |           (FEFamily Enum, cast to int/s)                    |
     |                                                             |
     |   end variable loop                                         |
     |                                                             |
     | 8.) The number of additional vectors (unsigned int),        |
     |                                                             |
     |    for each additional vector in the equation system object |
     |                                                             |
     |    9.) the name of the additional vector  (string)          |
     +-------------------------------------------------------------+
    
     end system loop
    
    
       for each system, handled through System::read_{serialized,parallel}_data():
       
     +--------------------------------------------------------------+
     | 10.) The global solution vector, re-ordered to be node-major |
     |     (More on this later.)                                    |
     |                                                              |
     |    for each additional vector in the equation system object  |
     |                                                              |
     |    11.) The global additional vector, re-ordered to be       |
     |         node-major (More on this later.)                     |
     +--------------------------------------------------------------+
    
     end system loop
   

Note that the actual IO is handled through the Xdr class (to be renamed later?) which provides a uniform interface to both the XDR (eXternal Data Representation) interface and standard ASCII output. Thus this one section of code will read XDR or ASCII files with no changes.

Definition at line 123 of file equation_systems_io.C.

References _mesh, _systems, add_system(), Xdr::close(), Xdr::data(), MeshTools::Private::fix_broken_node_and_element_numbering(), get_mesh(), get_system(), MeshTools::Private::globally_renumber_nodes_and_elements(), init(), local_file_name(), mesh, libMesh::processor_id(), read(), READ_ADDITIONAL_DATA, READ_DATA, System::read_header(), READ_HEADER, READ_LEGACY_FORMAT, Xdr::reading(), TRY_READ_IFEMS, and update().

Referenced by read().

00126 {
00190    // Set booleans from the read_flags argument
00191    const bool read_header          = read_flags & EquationSystems::READ_HEADER;
00192    const bool read_data            = read_flags & EquationSystems::READ_DATA;
00193    const bool read_additional_data = read_flags & EquationSystems::READ_ADDITIONAL_DATA;
00194    const bool read_legacy_format   = read_flags & EquationSystems::READ_LEGACY_FORMAT;
00195    const bool try_read_ifems       = read_flags & EquationSystems::TRY_READ_IFEMS;
00196          bool read_parallel_files  = false;
00197          
00198   // This will unzip a file with .bz2 as the extension, otherwise it
00199   // simply returns the name if the file need not be unzipped.
00200   Xdr io ((libMesh::processor_id() == 0) ? name : "", mode);
00201   libmesh_assert (io.reading());
00202           
00203   {
00204     // 1.)
00205     // Read the version header.
00206     std::string version = "legacy";
00207     if (!read_legacy_format)
00208       {
00209         if (libMesh::processor_id() == 0) io.data(version);     
00210         Parallel::broadcast(version);
00211         
00212         // All processors have the version header, if it does not contain
00213         // "libMesh" something then it is a legacy file.
00214         if (!(version.find("libMesh") < version.size()))
00215           {
00216             io.close();
00217             
00218             // Recursively call this read() function but with the 
00219             // EquationSystems::READ_LEGACY_FORMAT bit set.
00220             this->read (name, mode, (read_flags | EquationSystems::READ_LEGACY_FORMAT));
00221             return;
00222           }
00223 
00224         read_parallel_files = (version.rfind(" parallel") < version.size());
00225 
00226         // If requested that we try to read infinite element information,
00227         // and the string " with infinite elements" is not in the version,
00228         // then tack it on.  This is for compatibility reading ifem
00229         // files written prior to 11/10/2008 - BSK
00230         if (try_read_ifems)
00231           if (!(version.rfind(" with infinite elements") < version.size()))
00232             version += " with infinite elements";
00233         
00234       }
00235     else
00236       libmesh_deprecated();
00237 
00238   START_LOG("read()","EquationSystems");
00239       
00240     // 2.)  
00241     // Read the number of equation systems
00242     unsigned int n_sys=0;
00243     if (libMesh::processor_id() == 0) io.data (n_sys);
00244     Parallel::broadcast(n_sys);
00245 
00246     for (unsigned int sys=0; sys<n_sys; sys++)
00247       {
00248         // 3.)
00249         // Read the name of the sys-th equation system
00250         std::string sys_name;      
00251         if (libMesh::processor_id() == 0) io.data (sys_name);
00252         Parallel::broadcast(sys_name);
00253         
00254         // 4.)
00255         // Read the type of the sys-th equation system
00256         std::string sys_type;   
00257         if (libMesh::processor_id() == 0) io.data (sys_type);
00258         Parallel::broadcast(sys_type);
00259         
00260         if (read_header)
00261           this->add_system (sys_type, sys_name);
00262 
00263         // 5.) - 9.)
00264         // Let System::read_header() do the job
00265         System& new_system = this->get_system(sys_name);          
00266         new_system.read_header (io,
00267                                 version,
00268                                 read_header,
00269                                 read_additional_data,
00270                                 read_legacy_format);
00271       }
00272   }
00273       
00274 
00275 
00276   // Now we are ready to initialize the underlying data
00277   // structures. This will initialize the vectors for 
00278   // storage, the dof_map, etc...
00279   if (read_header) 
00280     this->init();
00281 
00282   // 10.) & 11.)
00283   // Read and set the numeric vector values
00284   if (read_data)
00285     {
00286       // the EquationSystems::read() method should look constant from the mesh
00287       // perspective, but we need to assign a temporary numbering to the nodes
00288       // and elements in the mesh, which requires that we abuse const_cast
00289       if (!read_legacy_format)
00290         {
00291           MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh());
00292           MeshTools::Private::globally_renumber_nodes_and_elements(mesh);
00293         }
00294 
00295       Xdr local_io (read_parallel_files ? local_file_name(name) : "", mode);
00296    
00297       std::map<std::string, System*>::iterator
00298         pos = _systems.begin();
00299       
00300       for (; pos != _systems.end(); ++pos)
00301         if (read_legacy_format)
00302           {
00303             libmesh_deprecated();
00304             pos->second->read_legacy_data (io, read_additional_data);
00305           }
00306         else
00307           if (read_parallel_files)
00308             pos->second->read_parallel_data   (local_io, read_additional_data);
00309           else
00310             pos->second->read_serialized_data (io, read_additional_data);
00311 
00312       
00313       // Undo the temporary numbering.
00314       if (!read_legacy_format)
00315         {
00316           if (dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh)))
00317             {
00318               ParallelMesh *mesh = dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh));    
00319               MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
00320             }
00321           else if (dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh)))
00322             {
00323               SerialMesh *mesh = dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh));    
00324               MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
00325             }
00326           else
00327             {
00328               std::cerr << "ERROR:  dynamic_cast<> to ParallelMesh and SerialMesh failed!"
00329                         << std::endl;
00330               libmesh_error();
00331             }     
00332         }
00333     }  
00334 
00335   STOP_LOG("read()","EquationSystems");
00336   
00337   // Localize each system's data
00338   this->update();
00339 }

template<typename T_sys >
T_sys & EquationSystems::add_system ( const std::string &  name  )  [inline]

Add the system named name to the systems array.

Definition at line 479 of file equation_systems.h.

References _add_system_to_nodes_and_elems(), _systems, and n_systems().

00480 {
00481   T_sys* ptr = NULL;
00482   
00483   if (!_systems.count(name))
00484     {
00485       ptr = new T_sys(*this, name, this->n_systems());
00486 
00487       _systems.insert (std::make_pair(name, ptr));   
00488 
00489       // Tell all the \p DofObject entities to add a system.
00490       this->_add_system_to_nodes_and_elems();
00491     }
00492   else
00493     {
00494       // We now allow redundant add_system calls, to make it
00495       // easier to load data from files for user-derived system
00496       // subclasses
00497 //      std::cerr << "ERROR: There was already a system"
00498 //              << " named " << name
00499 //              << std::endl;
00500 
00501 //      libmesh_error();
00502 
00503       ptr = &(this->get_system<T_sys>(name));
00504     }
00505 
00506   // Return a dynamically casted reference to the newly added System.
00507   return *ptr;
00508 }

System & EquationSystems::add_system ( const std::string &  system_type,
const std::string &  name 
)

Add the system of type system_type named name to the systems array.

Definition at line 304 of file equation_systems.C.

References get_system().

Referenced by _read_impl(), and ErrorVector::plot_error().

00306 {
00307   // Build a Newmark system
00308   if      (sys_type == "Newmark")
00309     this->add_system<NewmarkSystem> (name);
00310 
00311   // Build an Explicit system
00312   else if ((sys_type == "Explicit"))
00313     this->add_system<ExplicitSystem> (name);
00314 
00315   // Build an Implicit system
00316   else if ((sys_type == "Implicit") ||
00317            (sys_type == "Steady"  ))
00318     this->add_system<ImplicitSystem> (name);
00319 
00320   // build a transient implicit linear system
00321   else if ((sys_type == "Transient") ||
00322            (sys_type == "TransientImplicit") ||
00323            (sys_type == "TransientLinearImplicit"))
00324     this->add_system<TransientLinearImplicitSystem> (name);
00325 
00326   // build a transient implicit nonlinear system
00327   else if (sys_type == "TransientNonlinearImplicit")
00328     this->add_system<TransientNonlinearImplicitSystem> (name);
00329 
00330   // build a transient explicit system
00331   else if (sys_type == "TransientExplicit")
00332     this->add_system<TransientExplicitSystem> (name);
00333 
00334   // build a linear implicit system
00335   else if (sys_type == "LinearImplicit")
00336     this->add_system<LinearImplicitSystem> (name);
00337 
00338   // build a nonlinear implicit system
00339   else if (sys_type == "NonlinearImplicit")
00340     this->add_system<NonlinearImplicitSystem> (name);
00341 
00342 #if defined(LIBMESH_USE_COMPLEX_NUMBERS)
00343   // build a frequency system
00344   else if (sys_type == "Frequency")
00345     this->add_system<FrequencySystem> (name);
00346 #endif
00347 
00348   else
00349     {
00350       std::cerr << "ERROR: Unknown system type: "
00351                 << sys_type
00352                 << std::endl;
00353       libmesh_error();
00354     }
00355 
00356   // Return a reference to the new system
00357   //return (*this)(name);
00358   return this->get_system(name);
00359 }

void EquationSystems::adjoint_solve ( const QoISet qoi_indices = QoISet()  )  [virtual]

Call adjoint_solve on all the individual equation systems.

By default this function solves each system's adjoint once, in the reverse order from that in which they were added. For more sophisticated decoupled problems the user may with to override this behavior in a derived class.

Definition at line 405 of file equation_systems.C.

References get_system(), and n_systems().

00406 {
00407   libmesh_assert (this->n_systems());
00408 
00409   for (unsigned int i=this->n_systems(); i != 0; --i)
00410     this->get_system(i-1).adjoint_solve(qoi_indices);
00411 }

void EquationSystems::allgather (  ) 

Serializes a distributed mesh and its associated degree of freedom numbering for all systems

Definition at line 254 of file equation_systems.C.

References _mesh, MeshBase::allgather(), MeshBase::elements_begin(), MeshBase::elements_end(), get_system(), MeshBase::is_serial(), n_systems(), MeshBase::nodes_begin(), and MeshBase::nodes_end().

00255 {
00256   // A serial mesh means nothing needs to be done
00257   if (_mesh.is_serial())
00258     return;
00259 
00260   const unsigned int n_sys = this->n_systems();
00261 
00262   libmesh_assert (n_sys != 0);
00263 
00264   // Gather the mesh
00265   _mesh.allgather();
00266   
00267   // Tell all the \p DofObject entities how many systems
00268   // there are.
00269   {
00270     MeshBase::node_iterator       node_it  = _mesh.nodes_begin();
00271     const MeshBase::node_iterator node_end = _mesh.nodes_end();
00272 
00273     for ( ; node_it != node_end; ++node_it)
00274       (*node_it)->set_n_systems(n_sys);
00275     
00276     MeshBase::element_iterator       elem_it  = _mesh.elements_begin();
00277     const MeshBase::element_iterator elem_end = _mesh.elements_end();
00278     
00279     for ( ; elem_it != elem_end; ++elem_it)
00280       (*elem_it)->set_n_systems(n_sys);
00281   }
00282 
00283   // And distribute each system's dofs
00284   for (unsigned int i=0; i != this->n_systems(); ++i)
00285     this->get_system(i).get_dof_map().distribute_dofs(_mesh);
00286 }

void EquationSystems::build_discontinuous_solution_vector ( std::vector< Number > &  soln  )  const

Fill the input vector soln with solution values. The entries will be in variable-major format (corresponding to the names from build_variable_names()).

Definition at line 656 of file equation_systems.C.

References _mesh, _systems, MeshBase::active_elements_begin(), MeshBase::active_elements_end(), dim, DofMap::dof_indices(), System::get_dof_map(), Elem::infinite(), MeshBase::mesh_dimension(), Elem::n_nodes(), n_systems(), System::n_vars(), n_vars(), FEInterface::nodal_soln(), MeshBase::processor_id(), System::update_global_solution(), and System::variable_type().

Referenced by GMVIO::write_discontinuous_gmv().

00657 {
00658   libmesh_assert (this->n_systems());
00659 
00660   const unsigned int dim = _mesh.mesh_dimension();
00661   const unsigned int nv  = this->n_vars();
00662   unsigned int tw=0;
00663 
00664   // get the total weight
00665   {
00666     MeshBase::element_iterator       it  = _mesh.active_elements_begin();
00667     const MeshBase::element_iterator end = _mesh.active_elements_end(); 
00668 
00669     for ( ; it != end; ++it)
00670       tw += (*it)->n_nodes();
00671   }
00672   
00673 
00674   // Only if we are on processor zero, allocate the storage
00675   // to hold (number_of_nodes)*(number_of_variables) entries.
00676   if (_mesh.processor_id() == 0)
00677     soln.resize(tw*nv);
00678 
00679   std::vector<Number> sys_soln; 
00680 
00681   
00682   unsigned int var_num=0;
00683 
00684   // For each system in this EquationSystems object,
00685   // update the global solution and if we are on processor 0,
00686   // loop over the elements and build the nodal solution
00687   // from the element solution.  Then insert this nodal solution
00688   // into the vector passed to build_solution_vector.
00689   const_system_iterator       pos = _systems.begin();
00690   const const_system_iterator end = _systems.end();
00691 
00692   for (; pos != end; ++pos)
00693     {  
00694       const System& system  = *(pos->second);
00695       const unsigned int nv_sys = system.n_vars();
00696       
00697       system.update_global_solution (sys_soln, 0);
00698       
00699       if (_mesh.processor_id() == 0)
00700         {
00701           std::vector<Number>       elem_soln;   // The finite element solution
00702           std::vector<Number>       nodal_soln;  // The FE solution interpolated to the nodes
00703           std::vector<unsigned int> dof_indices; // The DOF indices for the finite element 
00704               
00705           for (unsigned int var=0; var<nv_sys; var++)
00706             {
00707               const FEType& fe_type    = system.variable_type(var);
00708 
00709               MeshBase::element_iterator       it  = _mesh.active_elements_begin();
00710               const MeshBase::element_iterator end = _mesh.active_elements_end(); 
00711 
00712               unsigned int nn=0;
00713               
00714               for ( ; it != end; ++it)
00715                 {
00716                   const Elem* elem = *it;
00717                   system.get_dof_map().dof_indices (elem, dof_indices, var);
00718                   
00719                   elem_soln.resize(dof_indices.size());
00720                   
00721                   for (unsigned int i=0; i<dof_indices.size(); i++)
00722                     elem_soln[i] = sys_soln[dof_indices[i]];
00723                   
00724                   FEInterface::nodal_soln (dim,
00725                                            fe_type,
00726                                            elem,
00727                                            elem_soln,
00728                                            nodal_soln);
00729 
00730 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00731                   // infinite elements should be skipped...
00732                   if (!elem->infinite())
00733 #endif
00734                     { 
00735                       libmesh_assert (nodal_soln.size() == elem->n_nodes());
00736                   
00737                       for (unsigned int n=0; n<elem->n_nodes(); n++)
00738                         {
00739                           soln[nv*(nn++) + (var + var_num)] +=
00740                             nodal_soln[n];
00741                         }
00742                     }
00743                 }
00744             }    
00745         }
00746 
00747       var_num += nv_sys;
00748     }
00749 }

void EquationSystems::build_solution_vector ( std::vector< Number > &  soln  )  const

Fill the input vector soln with solution values. The entries will be in variable-major format (corresponding to the names from build_variable_names()).

Definition at line 520 of file equation_systems.C.

References _mesh, _systems, MeshBase::active_local_elements_begin(), MeshBase::active_local_elements_end(), System::Variable::active_on_subdomain(), dim, DofMap::dof_indices(), System::get_dof_map(), System::Variable::implicitly_active(), Elem::infinite(), std::max(), MeshBase::max_node_id(), MeshBase::mesh_dimension(), Elem::n_nodes(), MeshBase::n_nodes(), n_systems(), System::n_vars(), n_vars(), FEInterface::nodal_soln(), Elem::node(), System::update_global_solution(), System::variable(), System::variable_type(), and libMesh::zero.

00521 {
00522   // This function must be run on all processors at once
00523   parallel_only();
00524 
00525   libmesh_assert (this->n_systems());
00526 
00527   const unsigned int dim = _mesh.mesh_dimension();
00528   const unsigned int nn  = _mesh.n_nodes();
00529   const unsigned int nv  = this->n_vars();
00530   const unsigned short int one = 1;
00531 
00532   // We'd better have a contiguous node numbering
00533   libmesh_assert (nn == _mesh.max_node_id());
00534 
00535   // allocate storage to hold
00536   // (number_of_nodes)*(number_of_variables) entries.
00537   soln.resize(nn*nv);
00538   
00539   // Zero out the soln vector
00540   std::fill (soln.begin(), soln.end(), libMesh::zero);
00541 
00542   std::vector<Number>  sys_soln;
00543 
00544   // (Note that we use an unsigned short int here even though an
00545   // unsigned char would be more that sufficient.  The MPI 1.1
00546   // standard does not require that MPI_SUM, MPI_PROD etc... be
00547   // implemented for char data types. 12/23/2003 - BSK)  
00548   std::vector<unsigned short int> node_conn(nn), repeat_count(nn);
00549 
00550   // Get the number of elements that share each node.  We will
00551   // compute the average value at each node.  This is particularly
00552   // useful for plotting discontinuous data.
00553   MeshBase::element_iterator       e_it  = _mesh.active_local_elements_begin();
00554   const MeshBase::element_iterator e_end = _mesh.active_local_elements_end(); 
00555 
00556   for ( ; e_it != e_end; ++e_it)
00557     for (unsigned int n=0; n<(*e_it)->n_nodes(); n++)
00558       node_conn[(*e_it)->node(n)]++;
00559 
00560   // Gather the distributed node_conn arrays in the case of
00561   // multiple processors.
00562   Parallel::sum(node_conn);
00563 
00564   unsigned int var_num=0;
00565 
00566   // For each system in this EquationSystems object,
00567   // update the global solution and if we are on processor 0,
00568   // loop over the elements and build the nodal solution
00569   // from the element solution.  Then insert this nodal solution
00570   // into the vector passed to build_solution_vector.
00571   const_system_iterator       pos = _systems.begin();
00572   const const_system_iterator end = _systems.end();
00573   
00574   for (; pos != end; ++pos)
00575     {  
00576       const System& system  = *(pos->second);
00577       const unsigned int nv_sys = system.n_vars();
00578       
00579       system.update_global_solution (sys_soln);
00580       
00581       std::vector<Number>       elem_soln;   // The finite element solution
00582       std::vector<Number>       nodal_soln;  // The FE solution interpolated to the nodes
00583       std::vector<unsigned int> dof_indices; // The DOF indices for the finite element 
00584       
00585       for (unsigned int var=0; var<nv_sys; var++)
00586         {
00587           const FEType& fe_type                   = system.variable_type(var);
00588           const System::Variable &var_description = system.variable(var);
00589           const DofMap &dof_map                   = system.get_dof_map();
00590 
00591           std::fill (repeat_count.begin(), repeat_count.end(), 0);
00592 
00593           MeshBase::element_iterator       it  = _mesh.active_local_elements_begin();
00594           const MeshBase::element_iterator end = _mesh.active_local_elements_end(); 
00595 
00596           for ( ; it != end; ++it)
00597             if (var_description.active_on_subdomain((*it)->subdomain_id()))
00598               {
00599                 const Elem* elem = *it;      
00600 
00601                 dof_map.dof_indices (elem, dof_indices, var);
00602               
00603                 elem_soln.resize(dof_indices.size());
00604               
00605                 for (unsigned int i=0; i<dof_indices.size(); i++)
00606                   elem_soln[i] = sys_soln[dof_indices[i]];
00607                   
00608                 FEInterface::nodal_soln (dim,
00609                                          fe_type,
00610                                          elem,
00611                                          elem_soln,
00612                                          nodal_soln);
00613 
00614 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00615                 // infinite elements should be skipped...
00616                 if (!elem->infinite())
00617 #endif
00618                   { 
00619                     libmesh_assert (nodal_soln.size() == elem->n_nodes());
00620                   
00621                     for (unsigned int n=0; n<elem->n_nodes(); n++)
00622                       {
00623                         repeat_count[elem->node(n)]++;
00624                         soln[nv*(elem->node(n)) + (var + var_num)] += nodal_soln[n];
00625                       }
00626                   }
00627               } // end loop over elements        
00628 
00629           // when a variable is active everywhere the repeat_count
00630           // and node_conn arrays should be identical, so let's
00631           // use that information to avoid unnecessary communication
00632           if (var_description.implicitly_active())
00633             repeat_count = node_conn;
00634             
00635           else
00636             Parallel::sum (repeat_count);
00637 
00638           for (unsigned int n=0; n<nn; n++)
00639             soln[nv*n + (var + var_num)] /= 
00640               static_cast<Real>(std::max (repeat_count[n], one));
00641           
00642       
00643         } // end loop on variables in this system
00644   
00645       var_num += nv_sys;
00646     } // end loop over systems
00647 
00648   // Now each processor has computed contriburions to the
00649   // soln vector.  Gather them all up.
00650   Parallel::sum(soln);
00651 }

void EquationSystems::build_solution_vector ( std::vector< Number > &  soln,
const std::string &  system_name,
const std::string &  variable_name = "all_vars" 
) const

Fill the input vector soln with the solution values for the system named name. Note that the input vector soln will only be assembled on processor 0, so this method is only applicable to outputting plot files from processor 0.

Definition at line 433 of file equation_systems.C.

Referenced by MeshOutput< MT >::_build_variable_names_and_solution_vector().

00436 {
00437   //TODO:[BSK] re-implement this from the method below
00438   libmesh_error();
00439 
00440 //   // Get a reference to the named system
00441 //   const System& system = this->get_system(system_name);
00442 
00443 //   // Get the number associated with the variable_name we are passed
00444 //   const unsigned short int variable_num = system.variable_number(variable_name);
00445 
00446 //   // Get the dimension of the current mesh
00447 //   const unsigned int dim = _mesh.mesh_dimension();
00448 
00449 //   // If we're on processor 0, allocate enough memory to hold the solution.
00450 //   // Since we're only looking at one variable, there will be one solution value
00451 //   // for each node in the mesh.
00452 //   if (_mesh.processor_id() == 0)
00453 //     soln.resize(_mesh.n_nodes());
00454 
00455 //   // Vector to hold the global solution from all processors
00456 //   std::vector<Number> sys_soln;
00457   
00458 //   // Update the global solution from all processors
00459 //   system.update_global_solution (sys_soln, 0);
00460   
00461 //   // Temporary vector to store the solution on an individual element.
00462 //   std::vector<Number>       elem_soln;   
00463 
00464 //   // The FE solution interpolated to the nodes
00465 //   std::vector<Number>       nodal_soln;  
00466 
00467 //   // The DOF indices for the element
00468 //   std::vector<unsigned int> dof_indices; 
00469 
00470 //   // Determine the finite/infinite element type used in this system 
00471 //   const FEType& fe_type    = system.variable_type(variable_num);
00472 
00473 //   // Define iterators to iterate over all the elements of the mesh
00474 //   const_active_elem_iterator       it (_mesh.elements_begin());
00475 //   const const_active_elem_iterator end(_mesh.elements_end());
00476 
00477 //   // Loop over elements
00478 //   for ( ; it != end; ++it)
00479 //     {
00480 //       // Convenient shortcut to the element pointer
00481 //       const Elem* elem = *it;
00482 
00483 //       // Fill the dof_indices vector for this variable
00484 //       system.get_dof_map().dof_indices(elem,
00485 //                                     dof_indices,
00486 //                                     variable_num);
00487 
00488 //       // Resize the element solution vector to fit the
00489 //       // dof_indices for this element.
00490 //       elem_soln.resize(dof_indices.size());
00491 
00492 //       // Transfer the system solution to the element
00493 //       // solution by mapping it through the dof_indices vector.
00494 //       for (unsigned int i=0; i<dof_indices.size(); i++)
00495 //      elem_soln[i] = sys_soln[dof_indices[i]];
00496 
00497 //       // Using the FE interface, compute the nodal_soln
00498 //       // for the current elemnt type given the elem_soln
00499 //       FEInterface::nodal_soln (dim,
00500 //                             fe_type,
00501 //                             elem,
00502 //                             elem_soln,
00503 //                             nodal_soln);
00504 
00505 //       // Sanity check -- make sure that there are the same number
00506 //       // of entries in the nodal_soln as there are nodes in the
00507 //       // element!
00508 //       libmesh_assert (nodal_soln.size() == elem->n_nodes());
00509 
00510 //       // Copy the nodal solution over into the correct place in
00511 //       // the global soln vector which will be returned to the user.
00512 //       for (unsigned int n=0; n<elem->n_nodes(); n++)
00513 //      soln[elem->node(n)] = nodal_soln[n];
00514 //     }
00515 }

void EquationSystems::build_variable_names ( std::vector< std::string > &  var_names  )  const

Fill the input vector var_names with the names of the variables for each system.

Definition at line 415 of file equation_systems.C.

References _systems, n_systems(), and n_vars().

Referenced by MeshOutput< MT >::_build_variable_names_and_solution_vector(), and GMVIO::write_discontinuous_gmv().

00416 {
00417   libmesh_assert (this->n_systems());
00418   
00419   var_names.resize (this->n_vars());
00420 
00421   unsigned int var_num=0;
00422   
00423   const_system_iterator       pos = _systems.begin();
00424   const const_system_iterator end = _systems.end();
00425 
00426   for (; pos != end; ++pos)
00427     for (unsigned int vn=0; vn<pos->second->n_vars(); vn++)
00428       var_names[var_num++] = pos->second->variable_name(vn);       
00429 }

void EquationSystems::clear (  ) 

Returns tha data structure to a pristine state.

Definition at line 68 of file equation_systems.C.

References _systems, Parameters::clear(), and parameters.

Referenced by ~EquationSystems().

00069 {
00070   // Clear any additional parameters
00071   parameters.clear ();
00072 
00073   // clear the systems.  We must delete them
00074   // since we newed them!
00075   while (!_systems.empty())
00076     {
00077       system_iterator pos = _systems.begin();
00078       
00079       System *sys = pos->second;
00080       delete sys;
00081       sys = NULL;
00082       
00083       _systems.erase (pos);
00084     }
00085 }

bool EquationSystems::compare ( const EquationSystems other_es,
const Real  threshold,
const bool  verbose 
) const

Returns:
true when this equation system contains identical data, up to the given threshold. Delegates most of the comparisons to perform to the responsible systems

Definition at line 753 of file equation_systems.C.

References _systems, System::compare(), get_system(), and n_systems().

00756 {
00757   // safety check, whether we handle at least the same number
00758   // of systems
00759   std::vector<bool> os_result;
00760 
00761   if (this->n_systems() != other_es.n_systems())
00762     {
00763       if (verbose)
00764         {
00765           std::cout << "  Fatal difference. This system handles " 
00766                     << this->n_systems() << " systems," << std::endl
00767                     << "  while the other system handles "
00768                     << other_es.n_systems() 
00769                     << " systems." << std::endl
00770                     << "  Aborting comparison." << std::endl;
00771         }
00772       return false;
00773     }
00774   else
00775     {
00776       // start comparing each system      
00777       const_system_iterator       pos = _systems.begin();
00778       const const_system_iterator end = _systems.end();
00779 
00780       for (; pos != end; ++pos)
00781         {
00782           const std::string& sys_name = pos->first;
00783           const System&  system        = *(pos->second);
00784       
00785           // get the other system
00786           const System& other_system   = other_es.get_system (sys_name);
00787 
00788           os_result.push_back (system.compare (other_system, threshold, verbose));
00789 
00790         }
00791 
00792     }
00793   
00794 
00795   // sum up the results
00796   if (os_result.size()==0)
00797     return true;
00798   else
00799     {
00800       bool os_identical;
00801       unsigned int n = 0;
00802           do
00803             {
00804               os_identical = os_result[n];
00805               n++;
00806             }
00807           while (os_identical && n<os_result.size());
00808           return os_identical;
00809         }
00810 }

void EquationSystems::delete_system ( const std::string &  name  ) 

Remove the system named name from the systems array. This function is now deprecated - write the libmesh-devel mailing list if you need it reimplemented.

Definition at line 366 of file equation_systems.C.

References _systems.

00367 {
00368   libmesh_deprecated();
00369 
00370   if (!_systems.count(name))
00371     {
00372       std::cerr << "ERROR: no system named "
00373                 << name  << std::endl;
00374 
00375       libmesh_error();
00376     }
00377   
00378   delete _systems[name];
00379   
00380   _systems.erase (name);
00381 }

std::string EquationSystems::get_info (  )  const

Returns:
a string containing information about the systems, flags, and parameters.

Definition at line 814 of file equation_systems.C.

References _systems, and n_systems().

Referenced by print_info().

00815 {
00816   std::ostringstream out;
00817   
00818   out << " EquationSystems\n"
00819       << "  n_systems()=" << this->n_systems() << '\n';
00820 
00821   // Print the info for the individual systems
00822   const_system_iterator       pos = _systems.begin();
00823   const const_system_iterator end = _systems.end();
00824 
00825   for (; pos != end; ++pos)
00826     out << pos->second->get_info();
00827 
00828   
00829 //   // Possibly print the parameters  
00830 //   if (!this->parameters.empty())
00831 //     {  
00832 //       out << "  n_parameters()=" << this->n_parameters() << '\n';
00833 //       out << "   Parameters:\n";
00834       
00835 //       for (std::map<std::string, Real>::const_iterator
00836 //           param = _parameters.begin(); param != _parameters.end();
00837 //         ++param)
00838 //      out << "    "
00839 //          << "\""
00840 //          << param->first
00841 //          << "\""
00842 //          << "="
00843 //          << param->second
00844 //          << '\n';
00845 //     }
00846   
00847   return out.str();
00848 }

MeshBase & EquationSystems::get_mesh (  )  [inline]

Returns:
a reference to the mesh

Definition at line 440 of file equation_systems.h.

References _mesh.

00441 {
00442   return _mesh;
00443 }

const MeshBase & EquationSystems::get_mesh (  )  const [inline]

MeshData & EquationSystems::get_mesh_data (  )  [inline]

Returns:
a reference to the mesh_data

Definition at line 455 of file equation_systems.h.

References _mesh_data.

00456 {
00457   libmesh_assert (_mesh_data != NULL);
00458   return *_mesh_data;
00459 }

const MeshData & EquationSystems::get_mesh_data (  )  const [inline]

Returns:
a constant reference to the mesh_data

Definition at line 447 of file equation_systems.h.

References _mesh_data.

00448 {
00449   libmesh_assert (_mesh_data != NULL);
00450   return *_mesh_data;
00451 }

System& EquationSystems::get_system ( const unsigned int  num  ) 

Returns:
a writeable referene to the system number num.

const System& EquationSystems::get_system ( const unsigned int  num  )  const

Returns:
a constant reference to system number num.

System& EquationSystems::get_system ( const std::string &  name  ) 

Returns:
a writeable referene to the system named name.

const System& EquationSystems::get_system ( const std::string &  name  )  const

Returns:
a constant reference to the system named name.

template<typename T_sys >
System & EquationSystems::get_system ( const unsigned int  num  )  [inline]

Returns:
a writeable referene to the system number num. The template argument defines the return type. For example, const SteadySystem& sys = eq.get_system<SteadySystem> (0); is an example of how the method might be used

Definition at line 565 of file equation_systems.h.

References _systems, and n_systems().

00566 {
00567   libmesh_assert (num < this->n_systems());
00568   
00569   const_system_iterator       pos = _systems.begin();
00570   const const_system_iterator end = _systems.end();
00571 
00572   for (; pos != end; ++pos)
00573     if (pos->second->number() == num)
00574       break;
00575 
00576   // Check for errors
00577   if (pos == end)
00578     {
00579       std::cerr << "ERROR: no system number " << num << " found!"
00580                 << std::endl;
00581       libmesh_error();
00582     }
00583 
00584   // Attempt dynamic cast
00585   T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
00586 
00587   // Check for failure of dynamic cast
00588   if (ptr == NULL)
00589     {
00590       LIBMESH_THROW(libMesh::DynamicCastFailure());
00591     }
00592 
00593   return *ptr; 
00594 }

template<typename T_sys >
const System & EquationSystems::get_system ( const unsigned int  num  )  const [inline]

Returns:
a constant reference to system number num. The template argument defines the return type. For example, const SteadySystem& sys = eq.get_system<SteadySystem> (0); is an example of how the method might be used

Definition at line 525 of file equation_systems.h.

References _systems, and n_systems().

00526 {
00527   libmesh_assert (num < this->n_systems());
00528 
00529 
00530   const_system_iterator       pos = _systems.begin();
00531   const const_system_iterator end = _systems.end();
00532 
00533   for (; pos != end; ++pos)
00534     if (pos->second->number() == num)
00535       break;
00536 
00537   // Check for errors
00538   if (pos == end)
00539     {
00540       std::cerr << "ERROR: no system number " << num << " found!"
00541                 << std::endl;
00542       libmesh_error();
00543     }
00544 
00545   // Attempt dynamic cast
00546   T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
00547 
00548   // Check for failure of dynamic cast
00549   if (ptr == NULL)
00550     {
00551       std::cerr << "ERROR: cannot convert system "
00552                 << num << " to requested type!"
00553                 << std::endl;
00554       libmesh_error();
00555     }
00556   
00557   return *ptr;
00558 }

template<typename T_sys >
System & EquationSystems::get_system ( const std::string &  name  )  [inline]

Returns:
a writeable referene to the system named name. The template argument defines the return type. For example, const SteadySystem& sys = eq.get_system<SteadySystem> ("sys"); is an example of how the method might be used

Definition at line 634 of file equation_systems.h.

References _systems.

00635 {
00636   system_iterator pos = _systems.find(name);
00637 
00638   // Check for errors
00639   if (pos == _systems.end())
00640     {
00641       std::cerr << "ERROR: no system named " << name << " found!"
00642                 << std::endl;
00643       libmesh_error();
00644     }
00645 
00646   // Attempt dynamic cast
00647   T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
00648 
00649   // Check for failure of dynamic cast
00650   if (ptr == NULL)
00651     {
00652       std::cerr << "ERROR: cannot convert system \""
00653                 << name << "\" to requested type!"
00654                 << std::endl;
00655       libmesh_error();
00656     }
00657 
00658   return *ptr; 
00659 }

template<typename T_sys >
const System & EquationSystems::get_system ( const std::string &  name  )  const [inline]

Returns:
a constant reference to the system named name. The template argument defines the return type. For example, const SteadySystem& sys = eq.get_system<SteadySystem> ("sys"); is an example of how the method might be used

Definition at line 603 of file equation_systems.h.

References _systems.

Referenced by ExactSolution::_compute_error(), _read_impl(), EnsightIO::add_scalar(), add_system(), EnsightIO::add_vector(), adjoint_solve(), allgather(), compare(), GMVIO::copy_nodal_solution(), ErrorEstimator::estimate_errors(), ExactSolution::ExactSolution(), init(), reinit(), sensitivity_solve(), VTKIO::solution_to_vtk(), solve(), VTKIO::system_vectors_to_vtk(), update(), EnsightIO::write_scalar_ascii(), and EnsightIO::write_vector_ascii().

00604 {
00605   const_system_iterator pos = _systems.find(name);
00606 
00607   // Check for errors
00608   if (pos == _systems.end())
00609     {
00610       std::cerr << "ERROR: no system named \"" << name << "\" found!"
00611                 << std::endl;
00612       libmesh_error();
00613     }
00614 
00615   // Attempt dynamic cast
00616   T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
00617 
00618   // Check for failure of dynamic cast
00619   if (ptr == NULL)
00620     {
00621       LIBMESH_THROW(libMesh::DynamicCastFailure());
00622     }
00623 
00624   return *ptr; 
00625 }

bool EquationSystems::has_mesh_data (  )  const [inline]

Returns:
true when the _mesh_data pointer is not NULL. This is needed because get_mesh_data will fail if it is NULL

Definition at line 462 of file equation_systems.h.

References _mesh_data.

00463 {
00464   return (_mesh_data!=NULL);
00465 }

bool EquationSystems::has_system ( const std::string &  name  )  const [inline]

Returns:
true if the system named name exists within this EquationSystems object.

Definition at line 513 of file equation_systems.h.

References _systems.

Referenced by EnsightIO::add_scalar(), and EnsightIO::add_vector().

00514 {
00515   if (_systems.find(name) == _systems.end())
00516     return false;
00517   return true;
00518 }

void EquationSystems::init (  ) 

Initialize all the systems

Definition at line 89 of file equation_systems.C.

References _mesh, MeshBase::delete_remote_elements(), MeshBase::elements_begin(), MeshBase::elements_end(), get_system(), libMesh::n_processors(), n_systems(), MeshBase::nodes_begin(), and MeshBase::nodes_end().

Referenced by _read_impl(), Solver::init(), and ErrorVector::plot_error().

00090 {
00091   const unsigned int n_sys = this->n_systems();
00092 
00093   libmesh_assert (n_sys != 0);
00094 
00095   // Distribute the mesh if possible
00096   if (libMesh::n_processors() > 1)
00097     _mesh.delete_remote_elements();
00098   
00099   // Tell all the \p DofObject entities how many systems
00100   // there are.
00101   {
00102     MeshBase::node_iterator       node_it  = _mesh.nodes_begin();
00103     const MeshBase::node_iterator node_end = _mesh.nodes_end();
00104 
00105     for ( ; node_it != node_end; ++node_it)
00106       (*node_it)->set_n_systems(n_sys);
00107     
00108     MeshBase::element_iterator       elem_it  = _mesh.elements_begin();
00109     const MeshBase::element_iterator elem_end = _mesh.elements_end();
00110     
00111     for ( ; elem_it != elem_end; ++elem_it)
00112       (*elem_it)->set_n_systems(n_sys);
00113   }
00114 
00115   for (unsigned int i=0; i != this->n_systems(); ++i)
00116     this->get_system(i).init();
00117 }

unsigned int EquationSystems::n_active_dofs (  )  const

Returns the number of active degrees of freedom for the EquationSystems object.

Definition at line 899 of file equation_systems.C.

References _systems.

00900 {
00901   unsigned int tot=0;
00902 
00903   const_system_iterator       pos = _systems.begin();
00904   const const_system_iterator end = _systems.end();
00905 
00906   for (; pos != end; ++pos)
00907     tot += pos->second->n_active_dofs();
00908 
00909   return tot;      
00910 }

unsigned int EquationSystems::n_dofs (  )  const

Returns:
the total number of degrees of freedom in all systems.

Definition at line 883 of file equation_systems.C.

References _systems.

00884 {
00885   unsigned int tot=0;
00886 
00887   const_system_iterator       pos = _systems.begin();
00888   const const_system_iterator end = _systems.end();
00889 
00890   for (; pos != end; ++pos)
00891     tot += pos->second->n_dofs();
00892 
00893   return tot;      
00894 }

unsigned int EquationSystems::n_vars (  )  const

Returns:
the total number of variables in all systems.

Definition at line 868 of file equation_systems.C.

References _systems.

Referenced by build_discontinuous_solution_vector(), build_solution_vector(), and build_variable_names().

00869 {
00870   unsigned int tot=0;
00871   
00872   const_system_iterator       pos = _systems.begin();
00873   const const_system_iterator end = _systems.end();
00874 
00875   for (; pos != end; ++pos)
00876     tot += pos->second->n_vars();
00877 
00878   return tot;
00879 }

void EquationSystems::print_info ( std::ostream &  os = std::cout  )  const

Prints information about the equation systems.

Definition at line 852 of file equation_systems.C.

References get_info().

Referenced by operator<<().

00853 {
00854   os << this->get_info()
00855      << std::endl;
00856 }

void EquationSystems::read ( const std::string &  name,
const libMeshEnums::XdrMODE  mode,
const unsigned int  read_flags = (READ_HEADER | READ_DATA) 
)

Read & initialize the systems from disk using the XDR data format. This format allows for machine-independent binary output.

Set which sections of the file to read by bitwise OR'ing the EquationSystems::ReadFlags enumeration together. For example, to read all sections of the file, set read_flags to: (READ_HEADER | READ_DATA | READ_ADDITIONAL_DATA)

Note that the equation system can be defined without initializing the data vectors to any solution values. This can be done by omitting READ_DATA in the read_flags parameter.

Definition at line 68 of file equation_systems_io.C.

References _read_impl(), and TRY_READ_IFEMS.

Referenced by _read_impl().

00071 {
00072 #ifdef LIBMESH_ENABLE_EXCEPTIONS
00073 
00074   // If we have exceptions enabled we can be considerate and try
00075   // to read old restart files which contain infinite element
00076   // information but do not have the " with infinite elements"
00077   // string in the version information.
00078   
00079   // First try the read the user requested
00080   try
00081     {
00082       this->_read_impl (name, mode, read_flags);
00083     }
00084 
00085   // If that fails, try it again but explicitly request we look for infinite element info      
00086   catch (...)
00087     {
00088       std::cout << "\n*********************************************************************\n"
00089                 << "READING THE FILE \"" << name << "\" FAILED.\n"
00090                 << "It is possible this file contains infinite element information,\n"
00091                 << "but the version string does not contain \" with infinite elements\"\n"
00092                 << "Let's try this again, but looking for infinite element information...\n"
00093                 << "*********************************************************************\n"
00094                 << std::endl;
00095 
00096       try
00097         {
00098           this->_read_impl (name, mode, read_flags | EquationSystems::TRY_READ_IFEMS);
00099         }
00100 
00101       // If all that failed, we are out of ideas here...
00102       catch (...)
00103         {
00104           std::cout << "\n*********************************************************************\n"
00105                     << "Well, at least we tried!\n"
00106                     << "Good Luck!!\n"
00107                     << "*********************************************************************\n"
00108                     << std::endl;
00109           throw;
00110         }
00111     }
00112 
00113 #else
00114 
00115   // no exceptions - cross your fingers...
00116   this->_read_impl (name, mode, read_flags);
00117   
00118 #endif // #ifdef LIBMESH_ENABLE_EXCEPTIONS
00119 }

void EquationSystems::reinit (  ) 

Reinitialize all the systems

Definition at line 121 of file equation_systems.C.

References _mesh, MeshRefinement::coarsen_elements(), MeshBase::contract(), DofMap::create_dof_constraints(), DofMap::distribute_dofs(), MeshBase::elements_begin(), MeshBase::elements_end(), MeshRefinement::face_level_mismatch_limit(), System::get_dof_map(), get_mesh(), get_system(), n_systems(), System::n_vars(), MeshBase::nodes_begin(), MeshBase::nodes_end(), DofMap::prepare_send_list(), DofMap::process_constraints(), System::prolong_vectors(), MeshRefinement::refine_elements(), System::restrict_vectors(), and System::user_constrain().

00122 {
00123   libmesh_assert (this->n_systems() != 0);
00124   
00125 #ifdef DEBUG
00126   // Make sure all the \p DofObject entities know how many systems
00127   // there are.
00128   {
00129     // All the nodes
00130     MeshBase::node_iterator       node_it  = _mesh.nodes_begin();
00131     const MeshBase::node_iterator node_end = _mesh.nodes_end();
00132 
00133     for ( ; node_it != node_end; ++node_it)
00134       libmesh_assert((*node_it)->n_systems() == this->n_systems());
00135     
00136     // All the elements
00137     MeshBase::element_iterator       elem_it  = _mesh.elements_begin();
00138     const MeshBase::element_iterator elem_end = _mesh.elements_end();
00139     
00140     for ( ; elem_it != elem_end; ++elem_it)
00141       libmesh_assert((*elem_it)->n_systems() == this->n_systems());
00142   }
00143 #endif
00144 
00145   // Localize each system's vectors
00146   for (unsigned int i=0; i != this->n_systems(); ++i)
00147     this->get_system(i).re_update();
00148 
00149 #ifdef LIBMESH_ENABLE_AMR
00150 
00151   bool dof_constraints_created = false;
00152   bool mesh_changed = false;
00153 
00154   // FIXME: For backwards compatibility, assume
00155   // refine_and_coarsen_elements or refine_uniformly have already
00156   // been called
00157   {
00158     for (unsigned int i=0; i != this->n_systems(); ++i)
00159       {
00160         System &sys = this->get_system(i);
00161 
00162         // Don't do anything if the system doesn't have any variables in it
00163         if(!sys.n_vars())
00164           continue;
00165         
00166         sys.get_dof_map().distribute_dofs(_mesh);
00167 
00168         // Recreate any hanging node constraints
00169         sys.get_dof_map().create_dof_constraints(_mesh);
00170 
00171         // Apply any user-defined constraints
00172         sys.user_constrain();
00173 
00174         // Expand any recursive constraints
00175         sys.get_dof_map().process_constraints();
00176 
00177         // And clean up the send_list before we use it again
00178         sys.get_dof_map().prepare_send_list();
00179 
00180         sys.prolong_vectors();
00181       }
00182     mesh_changed = true;
00183     dof_constraints_created = true;
00184   }
00185   
00186   // FIXME: Where should the user set maintain_level_one now??
00187   // Don't override previous settings, for now
00188 
00189   MeshRefinement mesh_refine(_mesh);
00190 
00191   mesh_refine.face_level_mismatch_limit() = false;
00192 
00193   // Try to coarsen the mesh, then restrict each system's vectors
00194   // if necessary
00195   if (mesh_refine.coarsen_elements())
00196     {
00197       for (unsigned int i=0; i != this->n_systems(); ++i)
00198         {
00199           System &sys = this->get_system(i);
00200           if (!dof_constraints_created)
00201             {
00202               sys.get_dof_map().distribute_dofs(_mesh);
00203               sys.get_dof_map().create_dof_constraints(_mesh);
00204               sys.user_constrain();
00205               sys.get_dof_map().process_constraints();
00206               sys.get_dof_map().prepare_send_list();
00207 
00208             }
00209           sys.restrict_vectors();
00210         }
00211       mesh_changed = true;
00212       dof_constraints_created = true;
00213     }
00214 
00215   // Once vectors are all restricted, we can delete
00216   // children of coarsened elements
00217   if (mesh_changed)
00218     this->get_mesh().contract();
00219   
00220   // Try to refine the mesh, then prolong each system's vectors
00221   // if necessary
00222   if (mesh_refine.refine_elements())
00223     {
00224       for (unsigned int i=0; i != this->n_systems(); ++i)
00225         {
00226           System &sys = this->get_system(i);
00227           if (!dof_constraints_created)
00228             {
00229               sys.get_dof_map().distribute_dofs(_mesh);
00230               sys.get_dof_map().create_dof_constraints(_mesh);
00231               sys.user_constrain();
00232               sys.get_dof_map().process_constraints();
00233               sys.get_dof_map().prepare_send_list();
00234 
00235             }
00236           sys.prolong_vectors();
00237         }
00238       mesh_changed = true;
00239       dof_constraints_created = true;
00240     }
00241 
00242   // If the mesh has changed, systems will need to create new dof
00243   // constraints and update their global solution vectors
00244   if (mesh_changed)
00245     {
00246       for (unsigned int i=0; i != this->n_systems(); ++i)
00247         this->get_system(i).reinit();
00248     }
00249 #endif // #ifdef LIBMESH_ENABLE_AMR
00250 }

void EquationSystems::sensitivity_solve ( const ParameterVector parameters  )  [virtual]

Call sensitivity_solve on all the individual equation systems.

By default this function solves each sensitivity system once, in the order in which in which they were added. For more sophisticated decoupled problems the user may with to override this behavior in a derived class.

Definition at line 395 of file equation_systems.C.

References get_system(), and n_systems().

00396 {
00397   libmesh_assert (this->n_systems());
00398 
00399   for (unsigned int i=0; i != this->n_systems(); ++i)
00400     this->get_system(i).sensitivity_solve(parameters);
00401 }

void EquationSystems::solve (  )  [virtual]

Call solve on all the individual equation systems.

By default this function solves each equation system once, in the order they were added. For more sophisticated decoupled problems the user may with to override this behavior in a derived class.

Definition at line 385 of file equation_systems.C.

References get_system(), and n_systems().

Referenced by Solver::solve().

00386 {
00387   libmesh_assert (this->n_systems());
00388 
00389   for (unsigned int i=0; i != this->n_systems(); ++i)
00390     this->get_system(i).solve();
00391 }

void EquationSystems::update (  ) 

Updates local values for all the systems

Definition at line 291 of file equation_systems.C.

References get_system(), and n_systems().

Referenced by _read_impl().

00292 {
00293   START_LOG("update()","EquationSystems");
00294   
00295   // Localize each system's vectors
00296   for (unsigned int i=0; i != this->n_systems(); ++i)
00297     this->get_system(i).update();
00298   
00299   STOP_LOG("update()","EquationSystems");
00300 }

void EquationSystems::write ( const std::string &  name,
const libMeshEnums::XdrMODE  mode,
const unsigned int  write_flags = (WRITE_DATA) 
) const

Write the systems to disk using the XDR data format. This format allows for machine-independent binary output.

Set the writing properties using the EquationSystems::WriteFlags enumeration. Set which sections to write out by bitwise OR'ing the enumeration values. Write everything by setting write_flags to: (WRITE_DATA | WRITE_ADDITIONAL_DATA)

Note that the solution data can be omitted by calling this routine with WRITE_DATA omitted in the write_flags argument.

This program implements the output of an EquationSystems object. This warrants some documentation. The output file essentially consists of 11 sections:

     1.) The version header.
     2.) The number of individual equation systems (unsigned int)
     
       for each system
                                                          
        3.)  The name of the system (string)            
        4.)  The type of the system (string)            
   
        handled through System::read():
   
     +-------------------------------------------------------------+
     |  5.) The number of variables in the system (unsigned int)   |
     |                                                             |
     |   for each variable in the system                           |
     |                                                             |
     |    6.) The name of the variable (string)                    |
     |                                                             |
     |    7.) Combined in an FEType:                               |
     |         - The approximation order(s) of the variable (Order |
     |           Enum, cast to int/s)                              |
     |         - The finite element family/ies of the variable     |
     |           (FEFamily Enum, cast to int/s)                    |
     |                                                             |
     |   end variable loop                                         |
     |                                                             |
     | 8.) The number of additional vectors (unsigned int),        |
     |                                                             |
     |    for each additional vector in the equation system object |
     |                                                             |
     |    9.) the name of the additional vector  (string)          |
     +-------------------------------------------------------------+
   
    end system loop
   
   
    for each system, handled through System::write_{serialized,parallel}_data():
       
     +--------------------------------------------------------------+
     | 10.) The global solution vector, re-ordered to be node-major |
     |     (More on this later.)                                    |
     |                                                              |
     |    for each additional vector in the equation system object  |
     |                                                              |
     |    11.) The global additional vector, re-ordered to be       |
     |         node-major (More on this later.)                     |
     +--------------------------------------------------------------+

    end system loop
   

Note that the actual IO is handled through the Xdr class (to be renamed later?) which provides a uniform interface to both the XDR (eXternal Data Representation) interface and standard ASCII output. Thus this one section of code will write XDR or ASCII files with no changes.

Definition at line 343 of file equation_systems_io.C.

References _mesh, _systems, Xdr::data(), MeshTools::Private::fix_broken_node_and_element_numbering(), get_mesh(), MeshTools::Private::globally_renumber_nodes_and_elements(), local_file_name(), mesh, n_systems(), libMesh::processor_id(), WRITE_ADDITIONAL_DATA, WRITE_DATA, WRITE_PARALLEL_FILES, and Xdr::writing().

00346 {
00410   // the EquationSystems::write() method should look constant,
00411   // but we need to assign a temporary numbering to the nodes
00412   // and elements in the mesh, which requires that we abuse const_cast
00413   {
00414     MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh());
00415     MeshTools::Private::globally_renumber_nodes_and_elements(mesh);
00416   }
00417   
00418    // set booleans from write_flags argument
00419    const bool write_data            = write_flags & EquationSystems::WRITE_DATA;
00420    const bool write_parallel_files  = write_flags & EquationSystems::WRITE_PARALLEL_FILES;
00421    const bool write_additional_data = write_flags & EquationSystems::WRITE_ADDITIONAL_DATA;
00422 
00423   // New scope so that io will close before we try to zip the file
00424   {
00425     Xdr io((libMesh::processor_id()==0) ? name : "", mode);    
00426     libmesh_assert (io.writing());
00427     
00428     START_LOG("write()","EquationSystems");
00429 
00430     const unsigned int proc_id = libMesh::processor_id();
00431     unsigned int n_sys         = this->n_systems();
00432     
00433     std::map<std::string, System*>::const_iterator
00434       pos = _systems.begin();
00435     
00436     std::string comment;
00437     char buf[256];
00438     
00439     // Only write the header information
00440     // if we are processor 0.
00441     if (proc_id == 0) 
00442       {
00443         // 1.)
00444         // Write the version header
00445         std::string version = "libMesh-0.7.0+";
00446         if (write_parallel_files) version += " parallel";
00447         
00448 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00449         version += " with infinite elements";
00450 #endif
00451         io.data (version, "# File Format Identifier");
00452         
00453         // 2.)  
00454         // Write the number of equation systems
00455         io.data (n_sys, "# No. of Equation Systems");
00456         
00457         while (pos != _systems.end())
00458           {
00459             // 3.)
00460             // Write the name of the sys_num-th system
00461             {
00462               const unsigned int sys_num = pos->second->number();
00463               std::string sys_name       = pos->first;
00464               
00465               comment =  "# Name, System No. ";
00466               std::sprintf(buf, "%d", sys_num);
00467               comment += buf;
00468           
00469               io.data (sys_name, comment.c_str());
00470             }
00471           
00472             // 4.)
00473             // Write the type of system handled
00474             {
00475               const unsigned int sys_num = pos->second->number();
00476               std::string sys_type       = pos->second->system_type();
00477 
00478               comment =  "# Type, System No. ";
00479               std::sprintf(buf, "%d", sys_num);
00480               comment += buf;
00481           
00482               io.data (sys_type, comment.c_str());
00483             }
00484         
00485             // 5.) - 9.)
00486             // Let System::write_header() do the job
00487             pos->second->write_header (io, version, write_additional_data);
00488             
00489             ++pos;
00490           }
00491       }
00492 
00493     // Start from the first system, again,
00494     // to write vectors to disk, if wanted
00495     if (write_data)
00496       {
00497         // open a parallel buffer if warranted.
00498         Xdr local_io (write_parallel_files ? local_file_name(name) : "", mode);
00499         
00500         for (pos = _systems.begin(); pos != _systems.end(); ++pos) 
00501           {
00502             // 10.) + 11.)
00503             if (write_parallel_files)
00504               pos->second->write_parallel_data (local_io,write_additional_data);
00505             else
00506               pos->second->write_serialized_data (io,write_additional_data);
00507           }
00508       }
00509         
00510     STOP_LOG("write()","EquationSystems");
00511   }
00512  
00513   // the EquationSystems::write() method should look constant,
00514   // but we need to undo the temporary numbering of the nodes
00515   // and elements in the mesh, which requires that we abuse const_cast
00516   if (dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh)))
00517     {
00518       ParallelMesh *mesh = dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh));    
00519       MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
00520     }
00521   else if (dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh)))
00522     {
00523       SerialMesh *mesh = dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh));    
00524       MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
00525     }
00526   else
00527     {
00528       std::cerr << "ERROR:  dynamic_cast<> to ParallelMesh and SerialMesh failed!"
00529                 << std::endl;
00530       libmesh_error();
00531     }
00532 }


Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const EquationSystems es 
) [friend]

Same as above, but allows you to also use stream syntax.

Definition at line 860 of file equation_systems.C.

00861 {
00862   es.print_info(os);
00863   return os;
00864 }


Member Data Documentation

A pointer to the MeshData object you would like to use. Can be NULL.

Definition at line 391 of file equation_systems.h.

Referenced by get_mesh_data(), and has_mesh_data().


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

Site Created By: libMesh Developers
Last modified: November 25 2009 03:44:11.

Hosted By:
SourceForge.net Logo