EquationSystems Class Reference
#include <equation_systems.h>

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 () |
| virtual void | clear () |
| virtual void | init () |
| virtual 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) |
| virtual 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 (const QoISet &qoi_indices=QoISet()) |
| virtual void | sensitivity_solve (const ParameterVector ¶meters) |
| 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 | read (const std::string &name, 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 |
| void | write (const std::string &name, const unsigned int write_flags=(WRITE_DATA)) const |
| virtual bool | compare (const EquationSystems &other_es, const Real threshold, const bool verbose) const |
| virtual std::string | get_info () const |
| void | print_info (std::ostream &os=libMesh::out) 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 () |
Static Public Member Functions | |
| static std::string | get_info () |
| static void | print_info (std::ostream &out=std::cout) |
| static unsigned int | n_objects () |
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 |
| typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| MeshBase & | _mesh |
| MeshData * | _mesh_data |
| std::map< std::string, System * > | _systems |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
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.
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 419 of file equation_systems.h.
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 106 of file reference_counter.h.
typedef std::map<std::string, System*>::iterator EquationSystems::system_iterator [protected] |
Typedef for system iterators
Definition at line 414 of file equation_systems.h.
Member Enumeration Documentation
Define enumeration to set properties in EquationSystems::read()
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()
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
Constructor.
Definition at line 50 of file equation_systems.C.
References parameters, and Parameters::set().
00051 : 00052 _mesh (m), 00053 _mesh_data (mesh_data) 00054 { 00055 // Set default parameters 00056 this->parameters.set<Real> ("linear solver tolerance") = TOLERANCE * TOLERANCE; 00057 this->parameters.set<unsigned int>("linear solver maximum iterations") = 5000; }
| 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 921 of file equation_systems.C.
References _mesh, MeshBase::elements_begin(), MeshBase::elements_end(), MeshBase::nodes_begin(), and MeshBase::nodes_end().
Referenced by add_system().
00922 { 00923 // All the nodes 00924 MeshBase::node_iterator node_it = _mesh.nodes_begin(); 00925 const MeshBase::node_iterator node_end = _mesh.nodes_end(); 00926 00927 for ( ; node_it != node_end; ++node_it) 00928 (*node_it)->add_system(); 00929 00930 // All the elements 00931 MeshBase::element_iterator elem_it = _mesh.elements_begin(); 00932 const MeshBase::element_iterator elem_end = _mesh.elements_end(); 00933 00934 for ( ; elem_it != elem_end; ++elem_it) 00935 (*elem_it)->add_system(); 00936 }
| 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 134 of file equation_systems_io.C.
References _mesh, _systems, add_system(), Parallel::broadcast(), Xdr::close(), Xdr::data(), libMesh::err, 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().
00137 { 00201 // Set booleans from the read_flags argument 00202 const bool read_header = read_flags & EquationSystems::READ_HEADER; 00203 const bool read_data = read_flags & EquationSystems::READ_DATA; 00204 const bool read_additional_data = read_flags & EquationSystems::READ_ADDITIONAL_DATA; 00205 const bool read_legacy_format = read_flags & EquationSystems::READ_LEGACY_FORMAT; 00206 const bool try_read_ifems = read_flags & EquationSystems::TRY_READ_IFEMS; 00207 bool read_parallel_files = false; 00208 00209 // This will unzip a file with .bz2 as the extension, otherwise it 00210 // simply returns the name if the file need not be unzipped. 00211 Xdr io ((libMesh::processor_id() == 0) ? name : "", mode); 00212 libmesh_assert (io.reading()); 00213 00214 { 00215 // 1.) 00216 // Read the version header. 00217 std::string version = "legacy"; 00218 if (!read_legacy_format) 00219 { 00220 if (libMesh::processor_id() == 0) io.data(version); 00221 Parallel::broadcast(version); 00222 00223 // All processors have the version header, if it does not contain 00224 // "libMesh" something then it is a legacy file. 00225 if (!(version.find("libMesh") < version.size())) 00226 { 00227 io.close(); 00228 00229 // Recursively call this read() function but with the 00230 // EquationSystems::READ_LEGACY_FORMAT bit set. 00231 this->read (name, mode, (read_flags | EquationSystems::READ_LEGACY_FORMAT)); 00232 return; 00233 } 00234 00235 read_parallel_files = (version.rfind(" parallel") < version.size()); 00236 00237 // If requested that we try to read infinite element information, 00238 // and the string " with infinite elements" is not in the version, 00239 // then tack it on. This is for compatibility reading ifem 00240 // files written prior to 11/10/2008 - BSK 00241 if (try_read_ifems) 00242 if (!(version.rfind(" with infinite elements") < version.size())) 00243 version += " with infinite elements"; 00244 00245 } 00246 else 00247 libmesh_deprecated(); 00248 00249 START_LOG("read()","EquationSystems"); 00250 00251 // 2.) 00252 // Read the number of equation systems 00253 unsigned int n_sys=0; 00254 if (libMesh::processor_id() == 0) io.data (n_sys); 00255 Parallel::broadcast(n_sys); 00256 00257 for (unsigned int sys=0; sys<n_sys; sys++) 00258 { 00259 // 3.) 00260 // Read the name of the sys-th equation system 00261 std::string sys_name; 00262 if (libMesh::processor_id() == 0) io.data (sys_name); 00263 Parallel::broadcast(sys_name); 00264 00265 // 4.) 00266 // Read the type of the sys-th equation system 00267 std::string sys_type; 00268 if (libMesh::processor_id() == 0) io.data (sys_type); 00269 Parallel::broadcast(sys_type); 00270 00271 if (read_header) 00272 this->add_system (sys_type, sys_name); 00273 00274 // 5.) - 9.) 00275 // Let System::read_header() do the job 00276 System& new_system = this->get_system(sys_name); 00277 new_system.read_header (io, 00278 version, 00279 read_header, 00280 read_additional_data, 00281 read_legacy_format); 00282 } 00283 } 00284 00285 00286 00287 // Now we are ready to initialize the underlying data 00288 // structures. This will initialize the vectors for 00289 // storage, the dof_map, etc... 00290 if (read_header) 00291 this->init(); 00292 00293 // 10.) & 11.) 00294 // Read and set the numeric vector values 00295 if (read_data) 00296 { 00297 // the EquationSystems::read() method should look constant from the mesh 00298 // perspective, but we need to assign a temporary numbering to the nodes 00299 // and elements in the mesh, which requires that we abuse const_cast 00300 if (!read_legacy_format) 00301 { 00302 MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh()); 00303 MeshTools::Private::globally_renumber_nodes_and_elements(mesh); 00304 } 00305 00306 Xdr local_io (read_parallel_files ? local_file_name(name) : "", mode); 00307 00308 std::map<std::string, System*>::iterator 00309 pos = _systems.begin(); 00310 00311 for (; pos != _systems.end(); ++pos) 00312 if (read_legacy_format) 00313 { 00314 libmesh_deprecated(); 00315 pos->second->read_legacy_data (io, read_additional_data); 00316 } 00317 else 00318 if (read_parallel_files) 00319 pos->second->read_parallel_data (local_io, read_additional_data); 00320 else 00321 pos->second->read_serialized_data (io, read_additional_data); 00322 00323 00324 // Undo the temporary numbering. 00325 if (!read_legacy_format) 00326 { 00327 if (dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh))) 00328 { 00329 ParallelMesh *mesh = dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh)); 00330 MeshTools::Private::fix_broken_node_and_element_numbering(*mesh); 00331 } 00332 else if (dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh))) 00333 { 00334 SerialMesh *mesh = dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh)); 00335 MeshTools::Private::fix_broken_node_and_element_numbering(*mesh); 00336 } 00337 else 00338 { 00339 libMesh::err << "ERROR: dynamic_cast<> to ParallelMesh and SerialMesh failed!" 00340 << std::endl; 00341 libmesh_error(); 00342 } 00343 } 00344 } 00345 00346 STOP_LOG("read()","EquationSystems"); 00347 00348 // Localize each system's data 00349 this->update(); 00350 }
| T_sys & EquationSystems::add_system | ( | const std::string & | name | ) | [inline] |
Add the system named name to the systems array.
Definition at line 492 of file equation_systems.h.
References _add_system_to_nodes_and_elems(), _systems, and n_systems().
00493 { 00494 T_sys* ptr = NULL; 00495 00496 if (!_systems.count(name)) 00497 { 00498 ptr = new T_sys(*this, name, this->n_systems()); 00499 00500 _systems.insert (std::make_pair(name, ptr)); 00501 00502 // Tell all the \p DofObject entities to add a system. 00503 this->_add_system_to_nodes_and_elems(); 00504 } 00505 else 00506 { 00507 // We now allow redundant add_system calls, to make it 00508 // easier to load data from files for user-derived system 00509 // subclasses 00510 // libMesh::err << "ERROR: There was already a system" 00511 // << " named " << name 00512 // << std::endl; 00513 00514 // libmesh_error(); 00515 00516 ptr = &(this->get_system<T_sys>(name)); 00517 } 00518 00519 // Return a dynamically casted reference to the newly added System. 00520 return *ptr; 00521 }
| System & EquationSystems::add_system | ( | const std::string & | system_type, | |
| const std::string & | name | |||
| ) | [virtual] |
Add the system of type system_type named name to the systems array.
Definition at line 304 of file equation_systems.C.
References _systems, libMesh::err, and get_system().
Referenced by _read_impl(), and ErrorVector::plot_error().
00306 { 00307 // If the user already built a system with this name, we'll 00308 // trust them and we'll use it. That way they can pre-add 00309 // non-standard derived system classes, and if their restart file 00310 // has some non-standard sys_type we won't throw an error. 00311 if (_systems.count(name)) 00312 { 00313 return this->get_system(name); 00314 } 00315 // Build a Newmark system 00316 else if (sys_type == "Newmark") 00317 this->add_system<NewmarkSystem> (name); 00318 00319 // Build an Explicit system 00320 else if ((sys_type == "Explicit")) 00321 this->add_system<ExplicitSystem> (name); 00322 00323 // Build an Implicit system 00324 else if ((sys_type == "Implicit") || 00325 (sys_type == "Steady" )) 00326 this->add_system<ImplicitSystem> (name); 00327 00328 // build a transient implicit linear system 00329 else if ((sys_type == "Transient") || 00330 (sys_type == "TransientImplicit") || 00331 (sys_type == "TransientLinearImplicit")) 00332 this->add_system<TransientLinearImplicitSystem> (name); 00333 00334 // build a transient implicit nonlinear system 00335 else if (sys_type == "TransientNonlinearImplicit") 00336 this->add_system<TransientNonlinearImplicitSystem> (name); 00337 00338 // build a transient explicit system 00339 else if (sys_type == "TransientExplicit") 00340 this->add_system<TransientExplicitSystem> (name); 00341 00342 // build a linear implicit system 00343 else if (sys_type == "LinearImplicit") 00344 this->add_system<LinearImplicitSystem> (name); 00345 00346 // build a nonlinear implicit system 00347 else if (sys_type == "NonlinearImplicit") 00348 this->add_system<NonlinearImplicitSystem> (name); 00349 00350 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) 00351 // build a frequency system 00352 else if (sys_type == "Frequency") 00353 this->add_system<FrequencySystem> (name); 00354 #endif 00355 00356 else 00357 { 00358 libMesh::err << "ERROR: Unknown system type: " 00359 << sys_type 00360 << std::endl; 00361 libmesh_error(); 00362 } 00363 00364 // Return a reference to the new system 00365 //return (*this)(name); 00366 return this->get_system(name); 00367 }
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 413 of file equation_systems.C.
References get_system(), and n_systems().
Referenced by AdjointRefinementEstimator::estimate_error().
00414 { 00415 libmesh_assert (this->n_systems()); 00416 00417 for (unsigned int i=this->n_systems(); i != 0; --i) 00418 this->get_system(i-1).adjoint_solve(qoi_indices); 00419 }
| 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 664 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 ExodusII_IO::write_discontinuous_exodusII(), and GMVIO::write_discontinuous_gmv().
00665 { 00666 libmesh_assert (this->n_systems()); 00667 00668 const unsigned int dim = _mesh.mesh_dimension(); 00669 const unsigned int nv = this->n_vars(); 00670 unsigned int tw=0; 00671 00672 // get the total weight 00673 { 00674 MeshBase::element_iterator it = _mesh.active_elements_begin(); 00675 const MeshBase::element_iterator end = _mesh.active_elements_end(); 00676 00677 for ( ; it != end; ++it) 00678 tw += (*it)->n_nodes(); 00679 } 00680 00681 00682 // Only if we are on processor zero, allocate the storage 00683 // to hold (number_of_nodes)*(number_of_variables) entries. 00684 if (_mesh.processor_id() == 0) 00685 soln.resize(tw*nv); 00686 00687 std::vector<Number> sys_soln; 00688 00689 00690 unsigned int var_num=0; 00691 00692 // For each system in this EquationSystems object, 00693 // update the global solution and if we are on processor 0, 00694 // loop over the elements and build the nodal solution 00695 // from the element solution. Then insert this nodal solution 00696 // into the vector passed to build_solution_vector. 00697 const_system_iterator pos = _systems.begin(); 00698 const const_system_iterator end = _systems.end(); 00699 00700 for (; pos != end; ++pos) 00701 { 00702 const System& system = *(pos->second); 00703 const unsigned int nv_sys = system.n_vars(); 00704 00705 system.update_global_solution (sys_soln, 0); 00706 00707 if (_mesh.processor_id() == 0) 00708 { 00709 std::vector<Number> elem_soln; // The finite element solution 00710 std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes 00711 std::vector<unsigned int> dof_indices; // The DOF indices for the finite element 00712 00713 for (unsigned int var=0; var<nv_sys; var++) 00714 { 00715 const FEType& fe_type = system.variable_type(var); 00716 00717 MeshBase::element_iterator it = _mesh.active_elements_begin(); 00718 const MeshBase::element_iterator end = _mesh.active_elements_end(); 00719 00720 unsigned int nn=0; 00721 00722 for ( ; it != end; ++it) 00723 { 00724 const Elem* elem = *it; 00725 system.get_dof_map().dof_indices (elem, dof_indices, var); 00726 00727 elem_soln.resize(dof_indices.size()); 00728 00729 for (unsigned int i=0; i<dof_indices.size(); i++) 00730 elem_soln[i] = sys_soln[dof_indices[i]]; 00731 00732 FEInterface::nodal_soln (dim, 00733 fe_type, 00734 elem, 00735 elem_soln, 00736 nodal_soln); 00737 00738 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00739 // infinite elements should be skipped... 00740 if (!elem->infinite()) 00741 #endif 00742 { 00743 libmesh_assert (nodal_soln.size() == elem->n_nodes()); 00744 00745 for (unsigned int n=0; n<elem->n_nodes(); n++) 00746 { 00747 soln[nv*(nn++) + (var + var_num)] += 00748 nodal_soln[n]; 00749 } 00750 } 00751 } 00752 } 00753 } 00754 00755 var_num += nv_sys; 00756 } 00757 }
| 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 528 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(), Parallel::sum(), System::update_global_solution(), System::variable(), System::variable_type(), and libMesh::zero.
00529 { 00530 // This function must be run on all processors at once 00531 parallel_only(); 00532 00533 libmesh_assert (this->n_systems()); 00534 00535 const unsigned int dim = _mesh.mesh_dimension(); 00536 const unsigned int nn = _mesh.n_nodes(); 00537 const unsigned int nv = this->n_vars(); 00538 const unsigned short int one = 1; 00539 00540 // We'd better have a contiguous node numbering 00541 libmesh_assert (nn == _mesh.max_node_id()); 00542 00543 // allocate storage to hold 00544 // (number_of_nodes)*(number_of_variables) entries. 00545 soln.resize(nn*nv); 00546 00547 // Zero out the soln vector 00548 std::fill (soln.begin(), soln.end(), libMesh::zero); 00549 00550 std::vector<Number> sys_soln; 00551 00552 // (Note that we use an unsigned short int here even though an 00553 // unsigned char would be more that sufficient. The MPI 1.1 00554 // standard does not require that MPI_SUM, MPI_PROD etc... be 00555 // implemented for char data types. 12/23/2003 - BSK) 00556 std::vector<unsigned short int> node_conn(nn), repeat_count(nn); 00557 00558 // Get the number of elements that share each node. We will 00559 // compute the average value at each node. This is particularly 00560 // useful for plotting discontinuous data. 00561 MeshBase::element_iterator e_it = _mesh.active_local_elements_begin(); 00562 const MeshBase::element_iterator e_end = _mesh.active_local_elements_end(); 00563 00564 for ( ; e_it != e_end; ++e_it) 00565 for (unsigned int n=0; n<(*e_it)->n_nodes(); n++) 00566 node_conn[(*e_it)->node(n)]++; 00567 00568 // Gather the distributed node_conn arrays in the case of 00569 // multiple processors. 00570 Parallel::sum(node_conn); 00571 00572 unsigned int var_num=0; 00573 00574 // For each system in this EquationSystems object, 00575 // update the global solution and if we are on processor 0, 00576 // loop over the elements and build the nodal solution 00577 // from the element solution. Then insert this nodal solution 00578 // into the vector passed to build_solution_vector. 00579 const_system_iterator pos = _systems.begin(); 00580 const const_system_iterator end = _systems.end(); 00581 00582 for (; pos != end; ++pos) 00583 { 00584 const System& system = *(pos->second); 00585 const unsigned int nv_sys = system.n_vars(); 00586 00587 system.update_global_solution (sys_soln); 00588 00589 std::vector<Number> elem_soln; // The finite element solution 00590 std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes 00591 std::vector<unsigned int> dof_indices; // The DOF indices for the finite element 00592 00593 for (unsigned int var=0; var<nv_sys; var++) 00594 { 00595 const FEType& fe_type = system.variable_type(var); 00596 const System::Variable &var_description = system.variable(var); 00597 const DofMap &dof_map = system.get_dof_map(); 00598 00599 std::fill (repeat_count.begin(), repeat_count.end(), 0); 00600 00601 MeshBase::element_iterator it = _mesh.active_local_elements_begin(); 00602 const MeshBase::element_iterator end = _mesh.active_local_elements_end(); 00603 00604 for ( ; it != end; ++it) 00605 if (var_description.active_on_subdomain((*it)->subdomain_id())) 00606 { 00607 const Elem* elem = *it; 00608 00609 dof_map.dof_indices (elem, dof_indices, var); 00610 00611 elem_soln.resize(dof_indices.size()); 00612 00613 for (unsigned int i=0; i<dof_indices.size(); i++) 00614 elem_soln[i] = sys_soln[dof_indices[i]]; 00615 00616 FEInterface::nodal_soln (dim, 00617 fe_type, 00618 elem, 00619 elem_soln, 00620 nodal_soln); 00621 00622 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00623 // infinite elements should be skipped... 00624 if (!elem->infinite()) 00625 #endif 00626 { 00627 libmesh_assert (nodal_soln.size() == elem->n_nodes()); 00628 00629 for (unsigned int n=0; n<elem->n_nodes(); n++) 00630 { 00631 repeat_count[elem->node(n)]++; 00632 soln[nv*(elem->node(n)) + (var + var_num)] += nodal_soln[n]; 00633 } 00634 } 00635 } // end loop over elements 00636 00637 // when a variable is active everywhere the repeat_count 00638 // and node_conn arrays should be identical, so let's 00639 // use that information to avoid unnecessary communication 00640 if (var_description.implicitly_active()) 00641 repeat_count = node_conn; 00642 00643 else 00644 Parallel::sum (repeat_count); 00645 00646 for (unsigned int n=0; n<nn; n++) 00647 soln[nv*n + (var + var_num)] /= 00648 static_cast<Real>(std::max (repeat_count[n], one)); 00649 00650 00651 } // end loop on variables in this system 00652 00653 var_num += nv_sys; 00654 } // end loop over systems 00655 00656 // Now each processor has computed contriburions to the 00657 // soln vector. Gather them all up. 00658 Parallel::sum(soln); 00659 }
| 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 441 of file equation_systems.C.
Referenced by MeshOutput< MT >::_build_variable_names_and_solution_vector().
00444 { 00445 //TODO:[BSK] re-implement this from the method below 00446 libmesh_error(); 00447 00448 // // Get a reference to the named system 00449 // const System& system = this->get_system(system_name); 00450 00451 // // Get the number associated with the variable_name we are passed 00452 // const unsigned short int variable_num = system.variable_number(variable_name); 00453 00454 // // Get the dimension of the current mesh 00455 // const unsigned int dim = _mesh.mesh_dimension(); 00456 00457 // // If we're on processor 0, allocate enough memory to hold the solution. 00458 // // Since we're only looking at one variable, there will be one solution value 00459 // // for each node in the mesh. 00460 // if (_mesh.processor_id() == 0) 00461 // soln.resize(_mesh.n_nodes()); 00462 00463 // // Vector to hold the global solution from all processors 00464 // std::vector<Number> sys_soln; 00465 00466 // // Update the global solution from all processors 00467 // system.update_global_solution (sys_soln, 0); 00468 00469 // // Temporary vector to store the solution on an individual element. 00470 // std::vector<Number> elem_soln; 00471 00472 // // The FE solution interpolated to the nodes 00473 // std::vector<Number> nodal_soln; 00474 00475 // // The DOF indices for the element 00476 // std::vector<unsigned int> dof_indices; 00477 00478 // // Determine the finite/infinite element type used in this system 00479 // const FEType& fe_type = system.variable_type(variable_num); 00480 00481 // // Define iterators to iterate over all the elements of the mesh 00482 // const_active_elem_iterator it (_mesh.elements_begin()); 00483 // const const_active_elem_iterator end(_mesh.elements_end()); 00484 00485 // // Loop over elements 00486 // for ( ; it != end; ++it) 00487 // { 00488 // // Convenient shortcut to the element pointer 00489 // const Elem* elem = *it; 00490 00491 // // Fill the dof_indices vector for this variable 00492 // system.get_dof_map().dof_indices(elem, 00493 // dof_indices, 00494 // variable_num); 00495 00496 // // Resize the element solution vector to fit the 00497 // // dof_indices for this element. 00498 // elem_soln.resize(dof_indices.size()); 00499 00500 // // Transfer the system solution to the element 00501 // // solution by mapping it through the dof_indices vector. 00502 // for (unsigned int i=0; i<dof_indices.size(); i++) 00503 // elem_soln[i] = sys_soln[dof_indices[i]]; 00504 00505 // // Using the FE interface, compute the nodal_soln 00506 // // for the current elemnt type given the elem_soln 00507 // FEInterface::nodal_soln (dim, 00508 // fe_type, 00509 // elem, 00510 // elem_soln, 00511 // nodal_soln); 00512 00513 // // Sanity check -- make sure that there are the same number 00514 // // of entries in the nodal_soln as there are nodes in the 00515 // // element! 00516 // libmesh_assert (nodal_soln.size() == elem->n_nodes()); 00517 00518 // // Copy the nodal solution over into the correct place in 00519 // // the global soln vector which will be returned to the user. 00520 // for (unsigned int n=0; n<elem->n_nodes(); n++) 00521 // soln[elem->node(n)] = nodal_soln[n]; 00522 // } 00523 }
| 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 423 of file equation_systems.C.
References _systems, n_systems(), and n_vars().
Referenced by MeshOutput< MT >::_build_variable_names_and_solution_vector(), ExodusII_IO::write_discontinuous_exodusII(), and GMVIO::write_discontinuous_gmv().
00424 { 00425 libmesh_assert (this->n_systems()); 00426 00427 var_names.resize (this->n_vars()); 00428 00429 unsigned int var_num=0; 00430 00431 const_system_iterator pos = _systems.begin(); 00432 const const_system_iterator end = _systems.end(); 00433 00434 for (; pos != end; ++pos) 00435 for (unsigned int vn=0; vn<pos->second->n_vars(); vn++) 00436 var_names[var_num++] = pos->second->variable_name(vn); 00437 }
| void EquationSystems::clear | ( | ) | [virtual] |
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 [virtual] |
- Returns:
truewhen 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 761 of file equation_systems.C.
References _systems, System::compare(), get_system(), n_systems(), and libMesh::out.
00764 { 00765 // safety check, whether we handle at least the same number 00766 // of systems 00767 std::vector<bool> os_result; 00768 00769 if (this->n_systems() != other_es.n_systems()) 00770 { 00771 if (verbose) 00772 { 00773 libMesh::out << " Fatal difference. This system handles " 00774 << this->n_systems() << " systems," << std::endl 00775 << " while the other system handles " 00776 << other_es.n_systems() 00777 << " systems." << std::endl 00778 << " Aborting comparison." << std::endl; 00779 } 00780 return false; 00781 } 00782 else 00783 { 00784 // start comparing each system 00785 const_system_iterator pos = _systems.begin(); 00786 const const_system_iterator end = _systems.end(); 00787 00788 for (; pos != end; ++pos) 00789 { 00790 const std::string& sys_name = pos->first; 00791 const System& system = *(pos->second); 00792 00793 // get the other system 00794 const System& other_system = other_es.get_system (sys_name); 00795 00796 os_result.push_back (system.compare (other_system, threshold, verbose)); 00797 00798 } 00799 00800 } 00801 00802 00803 // sum up the results 00804 if (os_result.size()==0) 00805 return true; 00806 else 00807 { 00808 bool os_identical; 00809 unsigned int n = 0; 00810 do 00811 { 00812 os_identical = os_result[n]; 00813 n++; 00814 } 00815 while (os_identical && n<os_result.size()); 00816 return os_identical; 00817 } 00818 }
| 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 374 of file equation_systems.C.
References _systems, and libMesh::err.
| std::string ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 45 of file reference_counter.C.
References ReferenceCounter::_counts, Quality::name(), and libMesh::out.
Referenced by ReferenceCounter::print_info().
00046 { 00047 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00048 00049 std::ostringstream out; 00050 00051 out << '\n' 00052 << " ---------------------------------------------------------------------------- \n" 00053 << "| Reference count information |\n" 00054 << " ---------------------------------------------------------------------------- \n"; 00055 00056 for (Counts::iterator it = _counts.begin(); 00057 it != _counts.end(); ++it) 00058 { 00059 const std::string name(it->first); 00060 const unsigned int creations = it->second.first; 00061 const unsigned int destructions = it->second.second; 00062 00063 out << "| " << name << " reference count information:\n" 00064 << "| Creations: " << creations << '\n' 00065 << "| Destructions: " << destructions << '\n'; 00066 } 00067 00068 out << " ---------------------------------------------------------------------------- \n"; 00069 00070 return out.str(); 00071 00072 #else 00073 00074 return ""; 00075 00076 #endif 00077 }
| std::string EquationSystems::get_info | ( | ) | const [virtual] |
- Returns:
- a string containing information about the systems, flags, and parameters.
Definition at line 822 of file equation_systems.C.
References _systems, n_systems(), and libMesh::out.
Referenced by print_info().
00823 { 00824 std::ostringstream out; 00825 00826 out << " EquationSystems\n" 00827 << " n_systems()=" << this->n_systems() << '\n'; 00828 00829 // Print the info for the individual systems 00830 const_system_iterator pos = _systems.begin(); 00831 const const_system_iterator end = _systems.end(); 00832 00833 for (; pos != end; ++pos) 00834 out << pos->second->get_info(); 00835 00836 00837 // // Possibly print the parameters 00838 // if (!this->parameters.empty()) 00839 // { 00840 // out << " n_parameters()=" << this->n_parameters() << '\n'; 00841 // out << " Parameters:\n"; 00842 00843 // for (std::map<std::string, Real>::const_iterator 00844 // param = _parameters.begin(); param != _parameters.end(); 00845 // ++param) 00846 // out << " " 00847 // << "\"" 00848 // << param->first 00849 // << "\"" 00850 // << "=" 00851 // << param->second 00852 // << '\n'; 00853 // } 00854 00855 return out.str(); 00856 }
| MeshBase & EquationSystems::get_mesh | ( | ) | [inline] |
- Returns:
- a reference to the mesh
Definition at line 453 of file equation_systems.h.
References _mesh.
00454 { 00455 return _mesh; 00456 }
| const MeshBase & EquationSystems::get_mesh | ( | ) | const [inline] |
- Returns:
- a constant reference to the mesh
Definition at line 445 of file equation_systems.h.
References _mesh.
Referenced by _read_impl(), AdjointRefinementEstimator::estimate_error(), MeshFunction::gradient(), MeshFunction::hessian(), MeshFunction::init(), MeshFunction::operator()(), reinit(), VTKIO::solution_to_vtk(), VTKIO::system_vectors_to_vtk(), write(), and VTKIO::write_equation_systems().
00446 { 00447 return _mesh; 00448 }
| MeshData & EquationSystems::get_mesh_data | ( | ) | [inline] |
- Returns:
- a reference to the mesh_data
Definition at line 468 of file equation_systems.h.
References _mesh_data.
00469 { 00470 libmesh_assert (_mesh_data != NULL); 00471 return *_mesh_data; 00472 }
| const MeshData & EquationSystems::get_mesh_data | ( | ) | const [inline] |
- Returns:
- a constant reference to the mesh_data
Definition at line 460 of file equation_systems.h.
References _mesh_data.
00461 { 00462 libmesh_assert (_mesh_data != NULL); 00463 return *_mesh_data; 00464 }
| 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.
| 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 578 of file equation_systems.h.
References _systems, libMesh::err, and n_systems().
00579 { 00580 libmesh_assert (num < this->n_systems()); 00581 00582 const_system_iterator pos = _systems.begin(); 00583 const const_system_iterator end = _systems.end(); 00584 00585 for (; pos != end; ++pos) 00586 if (pos->second->number() == num) 00587 break; 00588 00589 // Check for errors 00590 if (pos == end) 00591 { 00592 libMesh::err << "ERROR: no system number " << num << " found!" 00593 << std::endl; 00594 libmesh_error(); 00595 } 00596 00597 // Attempt dynamic cast 00598 T_sys* ptr = dynamic_cast<T_sys*>(pos->second); 00599 00600 // Check for failure of dynamic cast 00601 if (ptr == NULL) 00602 { 00603 LIBMESH_THROW(libMesh::DynamicCastFailure()); 00604 } 00605 00606 return *ptr; 00607 }
| 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 538 of file equation_systems.h.
References _systems, libMesh::err, and n_systems().
00539 { 00540 libmesh_assert (num < this->n_systems()); 00541 00542 00543 const_system_iterator pos = _systems.begin(); 00544 const const_system_iterator end = _systems.end(); 00545 00546 for (; pos != end; ++pos) 00547 if (pos->second->number() == num) 00548 break; 00549 00550 // Check for errors 00551 if (pos == end) 00552 { 00553 libMesh::err << "ERROR: no system number " << num << " found!" 00554 << std::endl; 00555 libmesh_error(); 00556 } 00557 00558 // Attempt dynamic cast 00559 T_sys* ptr = dynamic_cast<T_sys*>(pos->second); 00560 00561 // Check for failure of dynamic cast 00562 if (ptr == NULL) 00563 { 00564 libMesh::err << "ERROR: cannot convert system " 00565 << num << " to requested type!" 00566 << std::endl; 00567 libmesh_error(); 00568 } 00569 00570 return *ptr; 00571 }
| 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 647 of file equation_systems.h.
References _systems, and libMesh::err.
00648 { 00649 system_iterator pos = _systems.find(name); 00650 00651 // Check for errors 00652 if (pos == _systems.end()) 00653 { 00654 libMesh::err << "ERROR: no system named " << name << " found!" 00655 << std::endl; 00656 libmesh_error(); 00657 } 00658 00659 // Attempt dynamic cast 00660 T_sys* ptr = dynamic_cast<T_sys*>(pos->second); 00661 00662 // Check for failure of dynamic cast 00663 if (ptr == NULL) 00664 { 00665 std::cerr << "ERROR: cannot convert system \"" 00666 << name << "\" to requested type!" 00667 << std::endl; 00668 libmesh_error(); 00669 } 00670 00671 return *ptr; 00672 }
| 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 616 of file equation_systems.h.
References _systems, and libMesh::err.
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().
00617 { 00618 const_system_iterator pos = _systems.find(name); 00619 00620 // Check for errors 00621 if (pos == _systems.end()) 00622 { 00623 libMesh::err << "ERROR: no system named \"" << name << "\" found!" 00624 << std::endl; 00625 libmesh_error(); 00626 } 00627 00628 // Attempt dynamic cast 00629 T_sys* ptr = dynamic_cast<T_sys*>(pos->second); 00630 00631 // Check for failure of dynamic cast 00632 if (ptr == NULL) 00633 { 00634 LIBMESH_THROW(libMesh::DynamicCastFailure()); 00635 } 00636 00637 return *ptr; 00638 }
| 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 475 of file equation_systems.h.
References _mesh_data.
00476 { 00477 return (_mesh_data!=NULL); 00478 }
| bool EquationSystems::has_system | ( | const std::string & | name | ) | const [inline] |
- Returns:
- true if the system named
nameexists within this EquationSystems object.
Definition at line 526 of file equation_systems.h.
References _systems.
Referenced by EnsightIO::add_scalar(), and EnsightIO::add_vector().
| void ReferenceCounter::increment_constructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 150 of file reference_counter.h.
References ReferenceCounter::_counts, and Threads::spin_mtx.
Referenced by ReferenceCountedObject< Value >::ReferenceCountedObject().
00151 { 00152 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00153 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00154 00155 p.first++; 00156 }
| void ReferenceCounter::increment_destructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 168 of file reference_counter.h.
References ReferenceCounter::_counts, and Threads::spin_mtx.
Referenced by ReferenceCountedObject< Value >::~ReferenceCountedObject().
00169 { 00170 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00171 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00172 00173 p.second++; 00174 }
| void EquationSystems::init | ( | ) | [virtual] |
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 907 of file equation_systems.C.
References _systems.
00908 { 00909 unsigned int tot=0; 00910 00911 const_system_iterator pos = _systems.begin(); 00912 const const_system_iterator end = _systems.end(); 00913 00914 for (; pos != end; ++pos) 00915 tot += pos->second->n_active_dofs(); 00916 00917 return tot; 00918 }
| unsigned int EquationSystems::n_dofs | ( | ) | const |
- Returns:
- the total number of degrees of freedom in all systems.
Definition at line 891 of file equation_systems.C.
References _systems.
00892 { 00893 unsigned int tot=0; 00894 00895 const_system_iterator pos = _systems.begin(); 00896 const const_system_iterator end = _systems.end(); 00897 00898 for (; pos != end; ++pos) 00899 tot += pos->second->n_dofs(); 00900 00901 return tot; 00902 }
| static unsigned int ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 77 of file reference_counter.h.
References ReferenceCounter::_n_objects.
00078 { return _n_objects; }
| unsigned int EquationSystems::n_systems | ( | ) | const [inline] |
- Returns:
- the number of equation systems.
Definition at line 482 of file equation_systems.h.
References _systems.
Referenced by add_system(), adjoint_solve(), allgather(), build_discontinuous_solution_vector(), build_solution_vector(), build_variable_names(), compare(), GMVIO::copy_nodal_solution(), AdjointRefinementEstimator::estimate_error(), ErrorEstimator::estimate_errors(), ExactSolution::ExactSolution(), get_info(), get_system(), init(), reinit(), sensitivity_solve(), VTKIO::solution_to_vtk(), solve(), VTKIO::system_vectors_to_vtk(), update(), write(), and VTKIO::write_equation_systems().
00483 { 00484 return _systems.size(); 00485 }
| unsigned int EquationSystems::n_vars | ( | ) | const |
- Returns:
- the total number of variables in all systems.
Definition at line 876 of file equation_systems.C.
References _systems.
Referenced by build_discontinuous_solution_vector(), build_solution_vector(), and build_variable_names().
00877 { 00878 unsigned int tot=0; 00879 00880 const_system_iterator pos = _systems.begin(); 00881 const const_system_iterator end = _systems.end(); 00882 00883 for (; pos != end; ++pos) 00884 tot += pos->second->n_vars(); 00885 00886 return tot; 00887 }
| void ReferenceCounter::print_info | ( | std::ostream & | out = std::cout |
) | [static, inherited] |
Prints the reference information, by default to std::cout.
Definition at line 86 of file reference_counter.C.
References ReferenceCounter::get_info().
00087 { 00088 00089 out_stream << ReferenceCounter::get_info(); 00090 00091 }
| void EquationSystems::print_info | ( | std::ostream & | os = libMesh::out |
) | const |
Prints information about the equation systems, by default to libMesh::out.
Definition at line 860 of file equation_systems.C.
References get_info().
Referenced by operator<<().
00861 { 00862 os << this->get_info() 00863 << std::endl; 00864 }
| void EquationSystems::read | ( | const std::string & | name, | |
| const unsigned int | read_flags = (READ_HEADER | READ_DATA) | |||
| ) |
Definition at line 68 of file equation_systems_io.C.
References libMeshEnums::DECODE, read(), and libMeshEnums::READ.
| 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.
If XdrMODE is omitted, it will be inferred as READ for filenames containing .xda or as DECODE for filenames containing .xdr
Definition at line 79 of file equation_systems_io.C.
References _read_impl(), libMesh::out, and TRY_READ_IFEMS.
Referenced by _read_impl(), and read().
00082 { 00083 #ifdef LIBMESH_ENABLE_EXCEPTIONS 00084 00085 // If we have exceptions enabled we can be considerate and try 00086 // to read old restart files which contain infinite element 00087 // information but do not have the " with infinite elements" 00088 // string in the version information. 00089 00090 // First try the read the user requested 00091 try 00092 { 00093 this->_read_impl (name, mode, read_flags); 00094 } 00095 00096 // If that fails, try it again but explicitly request we look for infinite element info 00097 catch (...) 00098 { 00099 libMesh::out << "\n*********************************************************************\n" 00100 << "READING THE FILE \"" << name << "\" FAILED.\n" 00101 << "It is possible this file contains infinite element information,\n" 00102 << "but the version string does not contain \" with infinite elements\"\n" 00103 << "Let's try this again, but looking for infinite element information...\n" 00104 << "*********************************************************************\n" 00105 << std::endl; 00106 00107 try 00108 { 00109 this->_read_impl (name, mode, read_flags | EquationSystems::TRY_READ_IFEMS); 00110 } 00111 00112 // If all that failed, we are out of ideas here... 00113 catch (...) 00114 { 00115 libMesh::out << "\n*********************************************************************\n" 00116 << "Well, at least we tried!\n" 00117 << "Good Luck!!\n" 00118 << "*********************************************************************\n" 00119 << std::endl; 00120 throw; 00121 } 00122 } 00123 00124 #else 00125 00126 // no exceptions - cross your fingers... 00127 this->_read_impl (name, mode, read_flags); 00128 00129 #endif // #ifdef LIBMESH_ENABLE_EXCEPTIONS 00130 }
| void EquationSystems::reinit | ( | ) | [virtual] |
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().
Referenced by AdjointRefinementEstimator::estimate_error().
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 403 of file equation_systems.C.
References get_system(), and n_systems().
00404 { 00405 libmesh_assert (this->n_systems()); 00406 00407 for (unsigned int i=0; i != this->n_systems(); ++i) 00408 this->get_system(i).sensitivity_solve(parameters); 00409 }
| 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 393 of file equation_systems.C.
References get_system(), and n_systems().
Referenced by AdjointRefinementEstimator::estimate_error(), and Solver::solve().
00394 { 00395 libmesh_assert (this->n_systems()); 00396 00397 for (unsigned int i=0; i != this->n_systems(); ++i) 00398 this->get_system(i).solve(); 00399 }
| 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 unsigned int | write_flags = (WRITE_DATA) | |||
| ) | const |
Definition at line 354 of file equation_systems_io.C.
References libMeshEnums::ENCODE, write(), and libMeshEnums::WRITE.
| 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.
If XdrMODE is omitted, it will be inferred as WRITE for filenames containing .xda or as ENCODE for filenames containing .xdr
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 365 of file equation_systems_io.C.
References _mesh, _systems, Xdr::data(), libMesh::err, 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().
Referenced by write().
00368 { 00432 // the EquationSystems::write() method should look constant, 00433 // but we need to assign a temporary numbering to the nodes 00434 // and elements in the mesh, which requires that we abuse const_cast 00435 { 00436 MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh()); 00437 MeshTools::Private::globally_renumber_nodes_and_elements(mesh); 00438 } 00439 00440 // set booleans from write_flags argument 00441 const bool write_data = write_flags & EquationSystems::WRITE_DATA; 00442 const bool write_parallel_files = write_flags & EquationSystems::WRITE_PARALLEL_FILES; 00443 const bool write_additional_data = write_flags & EquationSystems::WRITE_ADDITIONAL_DATA; 00444 00445 // New scope so that io will close before we try to zip the file 00446 { 00447 Xdr io((libMesh::processor_id()==0) ? name : "", mode); 00448 libmesh_assert (io.writing()); 00449 00450 START_LOG("write()","EquationSystems"); 00451 00452 const unsigned int proc_id = libMesh::processor_id(); 00453 unsigned int n_sys = this->n_systems(); 00454 00455 std::map<std::string, System*>::const_iterator 00456 pos = _systems.begin(); 00457 00458 std::string comment; 00459 char buf[256]; 00460 00461 // Only write the header information 00462 // if we are processor 0. 00463 if (proc_id == 0) 00464 { 00465 // 1.) 00466 // Write the version header 00467 std::string version = "libMesh-0.7.0+"; 00468 if (write_parallel_files) version += " parallel"; 00469 00470 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00471 version += " with infinite elements"; 00472 #endif 00473 io.data (version, "# File Format Identifier"); 00474 00475 // 2.) 00476 // Write the number of equation systems 00477 io.data (n_sys, "# No. of Equation Systems"); 00478 00479 while (pos != _systems.end()) 00480 { 00481 // 3.) 00482 // Write the name of the sys_num-th system 00483 { 00484 const unsigned int sys_num = pos->second->number(); 00485 std::string sys_name = pos->first; 00486 00487 comment = "# Name, System No. "; 00488 std::sprintf(buf, "%d", sys_num); 00489 comment += buf; 00490 00491 io.data (sys_name, comment.c_str()); 00492 } 00493 00494 // 4.) 00495 // Write the type of system handled 00496 { 00497 const unsigned int sys_num = pos->second->number(); 00498 std::string sys_type = pos->second->system_type(); 00499 00500 comment = "# Type, System No. "; 00501 std::sprintf(buf, "%d", sys_num); 00502 comment += buf; 00503 00504 io.data (sys_type, comment.c_str()); 00505 } 00506 00507 // 5.) - 9.) 00508 // Let System::write_header() do the job 00509 pos->second->write_header (io, version, write_additional_data); 00510 00511 ++pos; 00512 } 00513 } 00514 00515 // Start from the first system, again, 00516 // to write vectors to disk, if wanted 00517 if (write_data) 00518 { 00519 // open a parallel buffer if warranted. 00520 Xdr local_io (write_parallel_files ? local_file_name(name) : "", mode); 00521 00522 for (pos = _systems.begin(); pos != _systems.end(); ++pos) 00523 { 00524 // 10.) + 11.) 00525 if (write_parallel_files) 00526 pos->second->write_parallel_data (local_io,write_additional_data); 00527 else 00528 pos->second->write_serialized_data (io,write_additional_data); 00529 } 00530 } 00531 00532 STOP_LOG("write()","EquationSystems"); 00533 } 00534 00535 // the EquationSystems::write() method should look constant, 00536 // but we need to undo the temporary numbering of the nodes 00537 // and elements in the mesh, which requires that we abuse const_cast 00538 if (dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh))) 00539 { 00540 ParallelMesh *mesh = dynamic_cast<ParallelMesh*>(const_cast<MeshBase*>(&_mesh)); 00541 MeshTools::Private::fix_broken_node_and_element_numbering(*mesh); 00542 } 00543 else if (dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh))) 00544 { 00545 SerialMesh *mesh = dynamic_cast<SerialMesh*>(const_cast<MeshBase*>(&_mesh)); 00546 MeshTools::Private::fix_broken_node_and_element_numbering(*mesh); 00547 } 00548 else 00549 { 00550 libMesh::err << "ERROR: dynamic_cast<> to ParallelMesh and SerialMesh failed!" 00551 << std::endl; 00552 libmesh_error(); 00553 } 00554 }
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 868 of file equation_systems.C.
00869 { 00870 es.print_info(os); 00871 return os; 00872 }
Member Data Documentation
ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 111 of file reference_counter.h.
Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().
MeshBase& EquationSystems::_mesh [protected] |
The mesh data structure
Definition at line 398 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().
MeshData* EquationSystems::_mesh_data [protected] |
A pointer to the MeshData object you would like to use. Can be NULL.
Definition at line 404 of file equation_systems.h.
Referenced by get_mesh_data(), and has_mesh_data().
Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 124 of file reference_counter.h.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 119 of file reference_counter.h.
Referenced by ReferenceCounter::n_objects(), ReferenceCounter::ReferenceCounter(), and ReferenceCounter::~ReferenceCounter().
std::map<std::string, System*> EquationSystems::_systems [protected] |
Data structure holding the systems.
Definition at line 409 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().
Data structure holding arbitrary parameters.
Definition at line 389 of file equation_systems.h.
Referenced by ExactSolution::_compute_error(), 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().
The documentation for this class was generated from the following files: