#include <equation_systems.h>
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 }
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 System & get_system (const std::string &name) const
System & get_system (const std::string &name)
const System & get_system (const unsigned int num) const
System & get_system (const unsigned int num)
System & add_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 ()
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 MeshBase & get_mesh () const
MeshBase & get_mesh ()
bool has_mesh_data () const
const MeshData & get_mesh_data () const
MeshData & get_mesh_data ()
void allgather ()
typedef std::map< std::string, System * >::iterator system_iterator
typedef std::map< std::string, System * >::const_iterator const_system_iterator
MeshBase & _mesh
MeshData * _mesh_data
std::map< std::string, System * > _systems
void _read_impl (const std::string &name, const libMeshEnums::XdrMODE, const unsigned int read_flags)
void _add_system_to_nodes_and_elems ()
std::ostream & operator<< (std::ostream &os, const EquationSystems &es)
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:
Definition at line 65 of file equation_systems.h.
Definition at line 396 of file equation_systems.h.
Definition at line 391 of file equation_systems.h.
Enumerator:
Definition at line 72 of file equation_systems.h.
{ READ_HEADER = 1,
READ_DATA = 2,
READ_ADDITIONAL_DATA = 4,
READ_LEGACY_FORMAT = 8,
TRY_READ_IFEMS = 16 };
Enumerator:
Definition at line 81 of file equation_systems.h.
{ WRITE_DATA = 1,
WRITE_ADDITIONAL_DATA = 2,
WRITE_PARALLEL_FILES = 4 };
Definition at line 50 of file equation_systems.C.
References parameters, and Parameters::set().
:
_mesh (m),
_mesh_data (mesh_data)
{
// Set default parameters
this->parameters.set<Real> ('linear solver tolerance') = TOLERANCE * TOLERANCE;
this->parameters.set<unsigned int>('linear solver maximum iterations') = 5000;
}
Definition at line 61 of file equation_systems.C.
References clear().
{
this->clear ();
}
Definition at line 903 of file equation_systems.C.
References _mesh, MeshBase::elements_begin(), MeshBase::elements_end(), MeshBase::nodes_begin(), and MeshBase::nodes_end().
Referenced by add_system().
{
// All the nodes
MeshBase::node_iterator node_it = _mesh.nodes_begin();
const MeshBase::node_iterator node_end = _mesh.nodes_end();
for ( ; node_it != node_end; ++node_it)
(*node_it)->add_system();
// All the elements
MeshBase::element_iterator elem_it = _mesh.elements_begin();
const MeshBase::element_iterator elem_end = _mesh.elements_end();
for ( ; elem_it != elem_end; ++elem_it)
(*elem_it)->add_system();
}
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(), 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().
{
// Set booleans from the read_flags argument
const bool read_header = read_flags & EquationSystems::READ_HEADER;
const bool read_data = read_flags & EquationSystems::READ_DATA;
const bool read_additional_data = read_flags & EquationSystems::READ_ADDITIONAL_DATA;
const bool read_legacy_format = read_flags & EquationSystems::READ_LEGACY_FORMAT;
const bool try_read_ifems = read_flags & EquationSystems::TRY_READ_IFEMS;
bool read_parallel_files = false;
// This will unzip a file with .bz2 as the extension, otherwise it
// simply returns the name if the file need not be unzipped.
Xdr io ((libMesh::processor_id() == 0) ? name : '', mode);
libmesh_assert (io.reading());
{
// 1.)
// Read the version header.
std::string version = 'legacy';
if (!read_legacy_format)
{
if (libMesh::processor_id() == 0) io.data(version);
Parallel::broadcast(version);
// All processors have the version header, if it does not contain
// 'libMesh' something then it is a legacy file.
if (!(version.find('libMesh') < version.size()))
{
io.close();
// Recursively call this read() function but with the
// EquationSystems::READ_LEGACY_FORMAT bit set.
this->read (name, mode, (read_flags | EquationSystems::READ_LEGACY_FORMAT));
return;
}
read_parallel_files = (version.rfind(' parallel') < version.size());
// If requested that we try to read infinite element information,
// and the string ' with infinite elements' is not in the version,
// then tack it on. This is for compatibility reading ifem
// files written prior to 11/10/2008 - BSK
if (try_read_ifems)
if (!(version.rfind(' with infinite elements') < version.size()))
version += ' with infinite elements';
}
else
libmesh_deprecated();
START_LOG('read()','EquationSystems');
// 2.)
// Read the number of equation systems
unsigned int n_sys=0;
if (libMesh::processor_id() == 0) io.data (n_sys);
Parallel::broadcast(n_sys);
for (unsigned int sys=0; sys<n_sys; sys++)
{
// 3.)
// Read the name of the sys-th equation system
std::string sys_name;
if (libMesh::processor_id() == 0) io.data (sys_name);
Parallel::broadcast(sys_name);
// 4.)
// Read the type of the sys-th equation system
std::string sys_type;
if (libMesh::processor_id() == 0) io.data (sys_type);
Parallel::broadcast(sys_type);
if (read_header)
this->add_system (sys_type, sys_name);
// 5.) - 9.)
// Let System::read_header() do the job
System& new_system = this->get_system(sys_name);
new_system.read_header (io,
version,
read_header,
read_additional_data,
read_legacy_format);
}
}
// Now we are ready to initialize the underlying data
// structures. This will initialize the vectors for
// storage, the dof_map, etc...
if (read_header)
this->init();
// 10.) & 11.)
// Read and set the numeric vector values
if (read_data)
{
// the EquationSystems::read() method should look constant from the mesh
// perspective, but we need to assign a temporary numbering to the nodes
// and elements in the mesh, which requires that we abuse const_cast
if (!read_legacy_format)
{
MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh());
MeshTools::Private::globally_renumber_nodes_and_elements(mesh);
}
Xdr local_io (read_parallel_files ? local_file_name(name) : '', mode);
std::map<std::string, System*>::iterator
pos = _systems.begin();
for (; pos != _systems.end(); ++pos)
if (read_legacy_format)
{
libmesh_deprecated();
pos->second->read_legacy_data (io, read_additional_data);
}
else
if (read_parallel_files)
pos->second->read_parallel_data (local_io, read_additional_data);
else
pos->second->read_serialized_data (io, read_additional_data);
// Undo the temporary numbering.
if (!read_legacy_format)
{
if (dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh)))
{
ParallelMesh *mesh = dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh));
MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
}
else if (dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh)))
{
SerialMesh *mesh = dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh));
MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
}
else
{
std::cerr << 'ERROR: dynamic_cast<> to ParallelMesh and SerialMesh failed!'
<< std::endl;
libmesh_error();
}
}
}
STOP_LOG('read()','EquationSystems');
// Localize each system's data
this->update();
}
Definition at line 304 of file equation_systems.C.
References get_system(), and Quality::name().
Referenced by _read_impl(), and ErrorVector::plot_error().
{
// Build a Newmark system
if (sys_type == 'Newmark')
this->add_system<NewmarkSystem> (name);
// Build an Explicit system
else if ((sys_type == 'Explicit'))
this->add_system<ExplicitSystem> (name);
// Build an Implicit system
else if ((sys_type == 'Implicit') ||
(sys_type == 'Steady' ))
this->add_system<ImplicitSystem> (name);
// build a transient implicit linear system
else if ((sys_type == 'Transient') ||
(sys_type == 'TransientImplicit') ||
(sys_type == 'TransientLinearImplicit'))
this->add_system<TransientLinearImplicitSystem> (name);
// build a transient implicit nonlinear system
else if (sys_type == 'TransientNonlinearImplicit')
this->add_system<TransientNonlinearImplicitSystem> (name);
// build a transient explicit system
else if (sys_type == 'TransientExplicit')
this->add_system<TransientExplicitSystem> (name);
// build a linear implicit system
else if (sys_type == 'LinearImplicit')
this->add_system<LinearImplicitSystem> (name);
// build a nonlinear implicit system
else if (sys_type == 'NonlinearImplicit')
this->add_system<NonlinearImplicitSystem> (name);
#if defined(LIBMESH_USE_COMPLEX_NUMBERS)
// build a frequency system
else if (sys_type == 'Frequency')
this->add_system<FrequencySystem> (name);
#endif
else
{
std::cerr << 'ERROR: Unknown system type: '
<< sys_type
<< std::endl;
libmesh_error();
}
// Return a reference to the new system
//return (*this)(name);
return this->get_system(name);
}
Definition at line 469 of file equation_systems.h.
References _add_system_to_nodes_and_elems(), _systems, n_systems(), and Quality::name().
{
T_sys* ptr = NULL;
if (!_systems.count(name))
{
ptr = new T_sys(*this, name, this->n_systems());
_systems.insert (std::make_pair(name, ptr));
// Tell all the
DofObject entities to add a system.
this->_add_system_to_nodes_and_elems();
}
else
{
// We now allow redundant add_system calls, to make it
// easier to load data from files for user-derived system
// subclasses
// std::cerr << 'ERROR: There was already a system'
// << ' named ' << name
// << std::endl;
// libmesh_error();
ptr = &(this->get_system<T_sys>(name));
}
// Return a dynamically casted reference to the newly added System.
return *ptr;
}
By default this function solves each equation system 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 395 of file equation_systems.C.
References get_system(), and n_systems().
Referenced by UniformRefinementEstimator::_estimate_error().
{
libmesh_assert (this->n_systems());
for (unsigned int i=this->n_systems(); i != 0; --i)
this->get_system(i-1).adjoint_solve();
}
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().
{
// A serial mesh means nothing needs to be done
if (_mesh.is_serial())
return;
const unsigned int n_sys = this->n_systems();
libmesh_assert (n_sys != 0);
// Gather the mesh
_mesh.allgather();
// Tell all the
DofObject entities how many systems
// there are.
{
MeshBase::node_iterator node_it = _mesh.nodes_begin();
const MeshBase::node_iterator node_end = _mesh.nodes_end();
for ( ; node_it != node_end; ++node_it)
(*node_it)->set_n_systems(n_sys);
MeshBase::element_iterator elem_it = _mesh.elements_begin();
const MeshBase::element_iterator elem_end = _mesh.elements_end();
for ( ; elem_it != elem_end; ++elem_it)
(*elem_it)->set_n_systems(n_sys);
}
// And distribute each system's dofs
for (unsigned int i=0; i != this->n_systems(); ++i)
this->get_system(i).get_dof_map().distribute_dofs(_mesh);
}
Definition at line 646 of file equation_systems.C.
References _mesh, _systems, MeshBase::active_elements_begin(), MeshBase::active_elements_end(), 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().
{
libmesh_assert (this->n_systems());
const unsigned int dim = _mesh.mesh_dimension();
const unsigned int nv = this->n_vars();
unsigned int tw=0;
// get the total weight
{
MeshBase::element_iterator it = _mesh.active_elements_begin();
const MeshBase::element_iterator end = _mesh.active_elements_end();
for ( ; it != end; ++it)
tw += (*it)->n_nodes();
}
// Only if we are on processor zero, allocate the storage
// to hold (number_of_nodes)*(number_of_variables) entries.
if (_mesh.processor_id() == 0)
soln.resize(tw*nv);
std::vector<Number> sys_soln;
unsigned int var_num=0;
// For each system in this EquationSystems object,
// update the global solution and if we are on processor 0,
// loop over the elements and build the nodal solution
// from the element solution. Then insert this nodal solution
// into the vector passed to build_solution_vector.
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
{
const System& system = *(pos->second);
const unsigned int nv_sys = system.n_vars();
system.update_global_solution (sys_soln, 0);
if (_mesh.processor_id() == 0)
{
std::vector<Number> elem_soln; // The finite element solution
std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes
std::vector<unsigned int> dof_indices; // The DOF indices for the finite element
for (unsigned int var=0; var<nv_sys; var++)
{
const FEType& fe_type = system.variable_type(var);
MeshBase::element_iterator it = _mesh.active_elements_begin();
const MeshBase::element_iterator end = _mesh.active_elements_end();
unsigned int nn=0;
for ( ; it != end; ++it)
{
const Elem* elem = *it;
system.get_dof_map().dof_indices (elem, dof_indices, var);
elem_soln.resize(dof_indices.size());
for (unsigned int i=0; i<dof_indices.size(); i++)
elem_soln[i] = sys_soln[dof_indices[i]];
FEInterface::nodal_soln (dim,
fe_type,
elem,
elem_soln,
nodal_soln);
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
// infinite elements should be skipped...
if (!elem->infinite())
#endif
{
libmesh_assert (nodal_soln.size() == elem->n_nodes());
for (unsigned int n=0; n<elem->n_nodes(); n++)
{
soln[nv*(nn++) + (var + var_num)] +=
nodal_soln[n];
}
}
}
}
}
var_num += nv_sys;
}
}
Definition at line 423 of file equation_systems.C.
Referenced by MeshOutput< MT >::_build_variable_names_and_solution_vector().
{
//TODO:[BSK] re-implement this from the method below
libmesh_error();
// // Get a reference to the named system
// const System& system = this->get_system(system_name);
// // Get the number associated with the variable_name we are passed
// const unsigned short int variable_num = system.variable_number(variable_name);
// // Get the dimension of the current mesh
// const unsigned int dim = _mesh.mesh_dimension();
// // If we're on processor 0, allocate enough memory to hold the solution.
// // Since we're only looking at one variable, there will be one solution value
// // for each node in the mesh.
// if (_mesh.processor_id() == 0)
// soln.resize(_mesh.n_nodes());
// // Vector to hold the global solution from all processors
// std::vector<Number> sys_soln;
// // Update the global solution from all processors
// system.update_global_solution (sys_soln, 0);
// // Temporary vector to store the solution on an individual element.
// std::vector<Number> elem_soln;
// // The FE solution interpolated to the nodes
// std::vector<Number> nodal_soln;
// // The DOF indices for the element
// std::vector<unsigned int> dof_indices;
// // Determine the finite/infinite element type used in this system
// const FEType& fe_type = system.variable_type(variable_num);
// // Define iterators to iterate over all the elements of the mesh
// const_active_elem_iterator it (_mesh.elements_begin());
// const const_active_elem_iterator end(_mesh.elements_end());
// // Loop over elements
// for ( ; it != end; ++it)
// {
// // Convenient shortcut to the element pointer
// const Elem* elem = *it;
// // Fill the dof_indices vector for this variable
// system.get_dof_map().dof_indices(elem,
// dof_indices,
// variable_num);
// // Resize the element solution vector to fit the
// // dof_indices for this element.
// elem_soln.resize(dof_indices.size());
// // Transfer the system solution to the element
// // solution by mapping it through the dof_indices vector.
// for (unsigned int i=0; i<dof_indices.size(); i++)
// elem_soln[i] = sys_soln[dof_indices[i]];
// // Using the FE interface, compute the nodal_soln
// // for the current elemnt type given the elem_soln
// FEInterface::nodal_soln (dim,
// fe_type,
// elem,
// elem_soln,
// nodal_soln);
// // Sanity check -- make sure that there are the same number
// // of entries in the nodal_soln as there are nodes in the
// // element!
// libmesh_assert (nodal_soln.size() == elem->n_nodes());
// // Copy the nodal solution over into the correct place in
// // the global soln vector which will be returned to the user.
// for (unsigned int n=0; n<elem->n_nodes(); n++)
// soln[elem->node(n)] = nodal_soln[n];
// }
}
Definition at line 510 of file equation_systems.C.
References _mesh, _systems, MeshBase::active_local_elements_begin(), MeshBase::active_local_elements_end(), System::Variable::active_on_subdomain(), 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.
{
// This function must be run on all processors at once
parallel_only();
libmesh_assert (this->n_systems());
const unsigned int dim = _mesh.mesh_dimension();
const unsigned int nn = _mesh.n_nodes();
const unsigned int nv = this->n_vars();
const unsigned short int one = 1;
// We'd better have a contiguous node numbering
libmesh_assert (nn == _mesh.max_node_id());
// allocate storage to hold
// (number_of_nodes)*(number_of_variables) entries.
soln.resize(nn*nv);
// Zero out the soln vector
std::fill (soln.begin(), soln.end(), libMesh::zero);
std::vector<Number> sys_soln;
// (Note that we use an unsigned short int here even though an
// unsigned char would be more that sufficient. The MPI 1.1
// standard does not require that MPI_SUM, MPI_PROD etc... be
// implemented for char data types. 12/23/2003 - BSK)
std::vector<unsigned short int> node_conn(nn), repeat_count(nn);
// Get the number of elements that share each node. We will
// compute the average value at each node. This is particularly
// useful for plotting discontinuous data.
MeshBase::element_iterator e_it = _mesh.active_local_elements_begin();
const MeshBase::element_iterator e_end = _mesh.active_local_elements_end();
for ( ; e_it != e_end; ++e_it)
for (unsigned int n=0; n<(*e_it)->n_nodes(); n++)
node_conn[(*e_it)->node(n)]++;
// Gather the distributed node_conn arrays in the case of
// multiple processors.
Parallel::sum(node_conn);
unsigned int var_num=0;
// For each system in this EquationSystems object,
// update the global solution and if we are on processor 0,
// loop over the elements and build the nodal solution
// from the element solution. Then insert this nodal solution
// into the vector passed to build_solution_vector.
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
{
const System& system = *(pos->second);
const unsigned int nv_sys = system.n_vars();
system.update_global_solution (sys_soln);
std::vector<Number> elem_soln; // The finite element solution
std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes
std::vector<unsigned int> dof_indices; // The DOF indices for the finite element
for (unsigned int var=0; var<nv_sys; var++)
{
const FEType& fe_type = system.variable_type(var);
const System::Variable &var_description = system.variable(var);
const DofMap &dof_map = system.get_dof_map();
std::fill (repeat_count.begin(), repeat_count.end(), 0);
MeshBase::element_iterator it = _mesh.active_local_elements_begin();
const MeshBase::element_iterator end = _mesh.active_local_elements_end();
for ( ; it != end; ++it)
if (var_description.active_on_subdomain((*it)->subdomain_id()))
{
const Elem* elem = *it;
dof_map.dof_indices (elem, dof_indices, var);
elem_soln.resize(dof_indices.size());
for (unsigned int i=0; i<dof_indices.size(); i++)
elem_soln[i] = sys_soln[dof_indices[i]];
FEInterface::nodal_soln (dim,
fe_type,
elem,
elem_soln,
nodal_soln);
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
// infinite elements should be skipped...
if (!elem->infinite())
#endif
{
libmesh_assert (nodal_soln.size() == elem->n_nodes());
for (unsigned int n=0; n<elem->n_nodes(); n++)
{
repeat_count[elem->node(n)]++;
soln[nv*(elem->node(n)) + (var + var_num)] += nodal_soln[n];
}
}
} // end loop over elements
// when a variable is active everywhere the repeat_count
// and node_conn arrays should be identical, so let's
// use that information to avoid unnecessary communication
if (var_description.implicitly_active())
repeat_count = node_conn;
else
Parallel::sum (repeat_count);
for (unsigned int n=0; n<nn; n++)
soln[nv*n + (var + var_num)] /=
static_cast<Real>(std::max (repeat_count[n], one));
} // end loop on variables in this system
var_num += nv_sys;
} // end loop over systems
// Now each processor has computed contriburions to the
// soln vector. Gather them all up.
Parallel::sum(soln);
}
Definition at line 405 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().
{
libmesh_assert (this->n_systems());
var_names.resize (this->n_vars());
unsigned int var_num=0;
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
for (unsigned int vn=0; vn<pos->second->n_vars(); vn++)
var_names[var_num++] = pos->second->variable_name(vn);
}
Definition at line 68 of file equation_systems.C.
References _systems, Parameters::clear(), and parameters.
Referenced by ~EquationSystems().
{
// Clear any additional parameters
parameters.clear ();
// clear the systems. We must delete them
// since we newed them!
while (!_systems.empty())
{
system_iterator pos = _systems.begin();
System *sys = pos->second;
delete sys;
sys = NULL;
_systems.erase (pos);
}
}
Definition at line 743 of file equation_systems.C.
References _systems, System::compare(), get_system(), and n_systems().
{
// safety check, whether we handle at least the same number
// of systems
std::vector<bool> os_result;
if (this->n_systems() != other_es.n_systems())
{
if (verbose)
{
std::cout << ' Fatal difference. This system handles '
<< this->n_systems() << ' systems,' << std::endl
<< ' while the other system handles '
<< other_es.n_systems()
<< ' systems.' << std::endl
<< ' Aborting comparison.' << std::endl;
}
return false;
}
else
{
// start comparing each system
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
{
const std::string& sys_name = pos->first;
const System& system = *(pos->second);
// get the other system
const System& other_system = other_es.get_system (sys_name);
os_result.push_back (system.compare (other_system, threshold, verbose));
}
}
// sum up the results
if (os_result.size()==0)
return true;
else
{
bool os_identical;
unsigned int n = 0;
do
{
os_identical = os_result[n];
n++;
}
while (os_identical && n<os_result.size());
return os_identical;
}
}
Definition at line 366 of file equation_systems.C.
References _systems, and Quality::name().
{
libmesh_deprecated();
if (!_systems.count(name))
{
std::cerr << 'ERROR: no system named '
<< name << std::endl;
libmesh_error();
}
delete _systems[name];
_systems.erase (name);
}
Definition at line 804 of file equation_systems.C.
References _systems, and n_systems().
Referenced by print_info().
{
std::ostringstream out;
out << ' EquationSystems
<< ' n_systems()=' << this->n_systems() << ';
// Print the info for the individual systems
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
out << pos->second->get_info();
// // Possibly print the parameters
// if (!this->parameters.empty())
// {
// out << ' n_parameters()=' << this->n_parameters() << ';
// out << ' Parameters:;
// for (std::map<std::string, Real>::const_iterator
// param = _parameters.begin(); param != _parameters.end();
// ++param)
// out << ' '
// << '''
// << param->first
// << '''
// << '='
// << param->second
// << ';
// }
return out.str();
}
Definition at line 422 of file equation_systems.h.
References _mesh.
Referenced by UniformRefinementEstimator::_estimate_error(), _read_impl(), MeshFunction::gradient(), MeshFunction::hessian(), MeshFunction::init(), MeshFunction::operator()(), reinit(), VTKIO::solution_to_vtk(), VTKIO::system_vectors_to_vtk(), write(), and VTKIO::write_equation_systems().
{
return _mesh;
}
Definition at line 430 of file equation_systems.h.
References _mesh.
{
return _mesh;
}
Definition at line 437 of file equation_systems.h.
References _mesh_data.
{
libmesh_assert (_mesh_data != NULL);
return *_mesh_data;
}
Definition at line 445 of file equation_systems.h.
References _mesh_data.
{
libmesh_assert (_mesh_data != NULL);
return *_mesh_data;
}
Definition at line 555 of file equation_systems.h.
References _systems, and n_systems().
{
libmesh_assert (num < this->n_systems());
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
if (pos->second->number() == num)
break;
// Check for errors
if (pos == end)
{
std::cerr << 'ERROR: no system number ' << num << ' found!'
<< std::endl;
libmesh_error();
}
// Attempt dynamic cast
T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
// Check for failure of dynamic cast
if (ptr == NULL)
{
LIBMESH_THROW(libMesh::DynamicCastFailure());
}
return *ptr;
}
Definition at line 515 of file equation_systems.h.
References _systems, and n_systems().
{
libmesh_assert (num < this->n_systems());
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
if (pos->second->number() == num)
break;
// Check for errors
if (pos == end)
{
std::cerr << 'ERROR: no system number ' << num << ' found!'
<< std::endl;
libmesh_error();
}
// Attempt dynamic cast
T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
// Check for failure of dynamic cast
if (ptr == NULL)
{
std::cerr << 'ERROR: cannot convert system '
<< num << ' to requested type!'
<< std::endl;
libmesh_error();
}
return *ptr;
}
Definition at line 624 of file equation_systems.h.
References _systems.
{
system_iterator pos = _systems.find(name);
// Check for errors
if (pos == _systems.end())
{
std::cerr << 'ERROR: no system named ' << name << ' found!'
<< std::endl;
libmesh_error();
}
// Attempt dynamic cast
T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
// Check for failure of dynamic cast
if (ptr == NULL)
{
std::cerr << 'ERROR: cannot convert system ''
<< name << '' to requested type!'
<< std::endl;
libmesh_error();
}
return *ptr;
}
Definition at line 593 of file equation_systems.h.
References _systems.
Referenced by ExactSolution::_compute_error(), UniformRefinementEstimator::_estimate_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(), VTKIO::solution_to_vtk(), solve(), VTKIO::system_vectors_to_vtk(), update(), EnsightIO::write_scalar_ascii(), and EnsightIO::write_vector_ascii().
{
const_system_iterator pos = _systems.find(name);
// Check for errors
if (pos == _systems.end())
{
std::cerr << 'ERROR: no system named '' << name << '' found!'
<< std::endl;
libmesh_error();
}
// Attempt dynamic cast
T_sys* ptr = dynamic_cast<T_sys*>(pos->second);
// Check for failure of dynamic cast
if (ptr == NULL)
{
LIBMESH_THROW(libMesh::DynamicCastFailure());
}
return *ptr;
}
Definition at line 452 of file equation_systems.h.
References _mesh_data.
{
return (_mesh_data!=NULL);
}
Definition at line 503 of file equation_systems.h.
References _systems.
Referenced by EnsightIO::add_scalar(), and EnsightIO::add_vector().
{
if (_systems.find(name) == _systems.end())
return false;
return true;
}
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().
{
const unsigned int n_sys = this->n_systems();
libmesh_assert (n_sys != 0);
// Distribute the mesh if possible
if (libMesh::n_processors() > 1)
_mesh.delete_remote_elements();
// Tell all the
DofObject entities how many systems
// there are.
{
MeshBase::node_iterator node_it = _mesh.nodes_begin();
const MeshBase::node_iterator node_end = _mesh.nodes_end();
for ( ; node_it != node_end; ++node_it)
(*node_it)->set_n_systems(n_sys);
MeshBase::element_iterator elem_it = _mesh.elements_begin();
const MeshBase::element_iterator elem_end = _mesh.elements_end();
for ( ; elem_it != elem_end; ++elem_it)
(*elem_it)->set_n_systems(n_sys);
}
for (unsigned int i=0; i != this->n_systems(); ++i)
this->get_system(i).init();
}
Definition at line 889 of file equation_systems.C.
References _systems.
{
unsigned int tot=0;
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
tot += pos->second->n_active_dofs();
return tot;
}
Definition at line 873 of file equation_systems.C.
References _systems.
{
unsigned int tot=0;
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
tot += pos->second->n_dofs();
return tot;
}
Definition at line 459 of file equation_systems.h.
References _systems.
Referenced by UniformRefinementEstimator::_estimate_error(), add_system(), adjoint_solve(), allgather(), build_discontinuous_solution_vector(), build_solution_vector(), build_variable_names(), compare(), GMVIO::copy_nodal_solution(), ErrorEstimator::estimate_errors(), ExactSolution::ExactSolution(), get_info(), get_system(), init(), reinit(), VTKIO::solution_to_vtk(), solve(), VTKIO::system_vectors_to_vtk(), update(), write(), and VTKIO::write_equation_systems().
{
return _systems.size();
}
Definition at line 858 of file equation_systems.C.
References _systems.
Referenced by build_discontinuous_solution_vector(), build_solution_vector(), and build_variable_names().
{
unsigned int tot=0;
const_system_iterator pos = _systems.begin();
const const_system_iterator end = _systems.end();
for (; pos != end; ++pos)
tot += pos->second->n_vars();
return tot;
}
Definition at line 842 of file equation_systems.C.
References get_info().
Referenced by operator<<().
{
os << this->get_info()
<< std::endl;
}
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().
{
#ifdef LIBMESH_ENABLE_EXCEPTIONS
// If we have exceptions enabled we can be considerate and try
// to read old restart files which contain infinite element
// information but do not have the ' with infinite elements'
// string in the version information.
// First try the read the user requested
try
{
this->_read_impl (name, mode, read_flags);
}
// If that fails, try it again but explicitly request we look for infinite element info
catch (...)
{
std::cout << '********************************************************************
<< 'READING THE FILE '' << name << '' FAILED.
<< 'It is possible this file contains infinite element information,
<< 'but the version string does not contain ' with infinite elements'
<< 'Let's try this again, but looking for infinite element information...
<< '*********************************************************************
<< std::endl;
try
{
this->_read_impl (name, mode, read_flags | EquationSystems::TRY_READ_IFEMS);
}
// If all that failed, we are out of ideas here...
catch (...)
{
std::cout << '********************************************************************
<< 'Well, at least we tried!
<< 'Good Luck!!
<< '*********************************************************************
<< std::endl;
throw;
}
}
#else
// no exceptions - cross your fingers...
this->_read_impl (name, mode, read_flags);
#endif // #ifdef LIBMESH_ENABLE_EXCEPTIONS
}
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().
Referenced by UniformRefinementEstimator::_estimate_error().
{
libmesh_assert (this->n_systems() != 0);
#ifdef DEBUG
// Make sure all the
DofObject entities know how many systems
// there are.
{
// All the nodes
MeshBase::node_iterator node_it = _mesh.nodes_begin();
const MeshBase::node_iterator node_end = _mesh.nodes_end();
for ( ; node_it != node_end; ++node_it)
libmesh_assert((*node_it)->n_systems() == this->n_systems());
// All the elements
MeshBase::element_iterator elem_it = _mesh.elements_begin();
const MeshBase::element_iterator elem_end = _mesh.elements_end();
for ( ; elem_it != elem_end; ++elem_it)
libmesh_assert((*elem_it)->n_systems() == this->n_systems());
}
#endif
// Localize each system's vectors
for (unsigned int i=0; i != this->n_systems(); ++i)
this->get_system(i).re_update();
#ifdef LIBMESH_ENABLE_AMR
bool dof_constraints_created = false;
bool mesh_changed = false;
// FIXME: For backwards compatibility, assume
// refine_and_coarsen_elements or refine_uniformly have already
// been called
{
for (unsigned int i=0; i != this->n_systems(); ++i)
{
System &sys = this->get_system(i);
// Don't do anything if the system doesn't have any variables in it
if(!sys.n_vars())
continue;
sys.get_dof_map().distribute_dofs(_mesh);
// Recreate any hanging node constraints
sys.get_dof_map().create_dof_constraints(_mesh);
// Apply any user-defined constraints
sys.user_constrain();
// Expand any recursive constraints
sys.get_dof_map().process_constraints();
// And clean up the send_list before we use it again
sys.get_dof_map().prepare_send_list();
sys.prolong_vectors();
}
mesh_changed = true;
dof_constraints_created = true;
}
// FIXME: Where should the user set maintain_level_one now??
// Don't override previous settings, for now
MeshRefinement mesh_refine(_mesh);
mesh_refine.face_level_mismatch_limit() = false;
// Try to coarsen the mesh, then restrict each system's vectors
// if necessary
if (mesh_refine.coarsen_elements())
{
for (unsigned int i=0; i != this->n_systems(); ++i)
{
System &sys = this->get_system(i);
if (!dof_constraints_created)
{
sys.get_dof_map().distribute_dofs(_mesh);
sys.get_dof_map().create_dof_constraints(_mesh);
sys.user_constrain();
sys.get_dof_map().process_constraints();
sys.get_dof_map().prepare_send_list();
}
sys.restrict_vectors();
}
mesh_changed = true;
dof_constraints_created = true;
}
// Once vectors are all restricted, we can delete
// children of coarsened elements
if (mesh_changed)
this->get_mesh().contract();
// Try to refine the mesh, then prolong each system's vectors
// if necessary
if (mesh_refine.refine_elements())
{
for (unsigned int i=0; i != this->n_systems(); ++i)
{
System &sys = this->get_system(i);
if (!dof_constraints_created)
{
sys.get_dof_map().distribute_dofs(_mesh);
sys.get_dof_map().create_dof_constraints(_mesh);
sys.user_constrain();
sys.get_dof_map().process_constraints();
sys.get_dof_map().prepare_send_list();
}
sys.prolong_vectors();
}
mesh_changed = true;
dof_constraints_created = true;
}
// If the mesh has changed, systems will need to create new dof
// constraints and update their global solution vectors
if (mesh_changed)
{
for (unsigned int i=0; i != this->n_systems(); ++i)
this->get_system(i).reinit();
}
#endif // #ifdef LIBMESH_ENABLE_AMR
}
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 UniformRefinementEstimator::_estimate_error(), and Solver::solve().
{
libmesh_assert (this->n_systems());
for (unsigned int i=0; i != this->n_systems(); ++i)
this->get_system(i).solve();
}
Definition at line 291 of file equation_systems.C.
References get_system(), and n_systems().
Referenced by _read_impl().
{
START_LOG('update()','EquationSystems');
// Localize each system's vectors
for (unsigned int i=0; i != this->n_systems(); ++i)
this->get_system(i).update();
STOP_LOG('update()','EquationSystems');
}
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(), n_systems(), libMesh::processor_id(), WRITE_ADDITIONAL_DATA, WRITE_DATA, WRITE_PARALLEL_FILES, and Xdr::writing().
{
// the EquationSystems::write() method should look constant,
// but we need to assign a temporary numbering to the nodes
// and elements in the mesh, which requires that we abuse const_cast
{
MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh());
MeshTools::Private::globally_renumber_nodes_and_elements(mesh);
}
// set booleans from write_flags argument
const bool write_data = write_flags & EquationSystems::WRITE_DATA;
const bool write_parallel_files = write_flags & EquationSystems::WRITE_PARALLEL_FILES;
const bool write_additional_data = write_flags & EquationSystems::WRITE_ADDITIONAL_DATA;
// New scope so that io will close before we try to zip the file
{
Xdr io((libMesh::processor_id()==0) ? name : '', mode);
libmesh_assert (io.writing());
START_LOG('write()','EquationSystems');
const unsigned int proc_id = libMesh::processor_id();
unsigned int n_sys = this->n_systems();
std::map<std::string, System*>::const_iterator
pos = _systems.begin();
std::string comment;
char buf[256];
// Only write the header information
// if we are processor 0.
if (proc_id == 0)
{
// 1.)
// Write the version header
std::string version = 'libMesh-0.7.0+';
if (write_parallel_files) version += ' parallel';
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
version += ' with infinite elements';
#endif
io.data (version, '# File Format Identifier');
// 2.)
// Write the number of equation systems
io.data (n_sys, '# No. of Equation Systems');
while (pos != _systems.end())
{
// 3.)
// Write the name of the sys_num-th system
{
const unsigned int sys_num = pos->second->number();
std::string sys_name = pos->first;
comment = '# Name, System No. ';
std::sprintf(buf, '%d', sys_num);
comment += buf;
io.data (sys_name, comment.c_str());
}
// 4.)
// Write the type of system handled
{
const unsigned int sys_num = pos->second->number();
std::string sys_type = pos->second->system_type();
comment = '# Type, System No. ';
std::sprintf(buf, '%d', sys_num);
comment += buf;
io.data (sys_type, comment.c_str());
}
// 5.) - 9.)
// Let System::write_header() do the job
pos->second->write_header (io, version, write_additional_data);
++pos;
}
}
// Start from the first system, again,
// to write vectors to disk, if wanted
if (write_data)
{
// open a parallel buffer if warranted.
Xdr local_io (write_parallel_files ? local_file_name(name) : '', mode);
for (pos = _systems.begin(); pos != _systems.end(); ++pos)
{
// 10.) + 11.)
if (write_parallel_files)
pos->second->write_parallel_data (local_io,write_additional_data);
else
pos->second->write_serialized_data (io,write_additional_data);
}
}
STOP_LOG('write()','EquationSystems');
}
// the EquationSystems::write() method should look constant,
// but we need to undo the temporary numbering of the nodes
// and elements in the mesh, which requires that we abuse const_cast
if (dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh)))
{
ParallelMesh *mesh = dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh));
MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
}
else if (dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh)))
{
SerialMesh *mesh = dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh));
MeshTools::Private::fix_broken_node_and_element_numbering(*mesh);
}
else
{
std::cerr << 'ERROR: dynamic_cast<> to ParallelMesh and SerialMesh failed!'
<< std::endl;
libmesh_error();
}
}
Definition at line 850 of file equation_systems.C.
{
es.print_info(os);
return os;
}
Definition at line 375 of file equation_systems.h.
Referenced by _add_system_to_nodes_and_elems(), _read_impl(), allgather(), build_discontinuous_solution_vector(), build_solution_vector(), get_mesh(), init(), reinit(), and write().
Definition at line 381 of file equation_systems.h.
Referenced by get_mesh_data(), and has_mesh_data().
Definition at line 386 of file equation_systems.h.
Referenced by _read_impl(), add_system(), build_discontinuous_solution_vector(), build_solution_vector(), build_variable_names(), clear(), compare(), delete_system(), get_info(), get_system(), has_system(), n_active_dofs(), n_dofs(), n_systems(), n_vars(), and write().
Definition at line 366 of file equation_systems.h.
Referenced by ExactSolution::_compute_error(), LinearImplicitSystem::adjoint_solve(), NewmarkSystem::clear(), clear(), FrequencySystem::clear_all(), EquationSystems(), ExactErrorEstimator::find_squared_element_error(), FEComputeData::init(), FrequencySystem::init_data(), FrequencySystem::n_frequencies(), NewmarkSystem::NewmarkSystem(), NonlinearImplicitSystem::NonlinearImplicitSystem(), FrequencySystem::set_current_frequency(), FrequencySystem::set_frequencies(), FrequencySystem::set_frequencies_by_range(), FrequencySystem::set_frequencies_by_steps(), NewmarkSystem::set_newmark_parameters(), NonlinearImplicitSystem::set_solver_parameters(), LinearImplicitSystem::solve(), FrequencySystem::solve(), and EigenSystem::solve().
Generated automatically by Doxygen for libMesh from the source code.