libMesh::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, READ_BASIC_ONLY = 32 } |
| enum | WriteFlags { WRITE_DATA = 1, WRITE_ADDITIONAL_DATA = 2, WRITE_PARALLEL_FILES = 4, WRITE_SERIAL_FILES = 8 } |
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 |
| std::size_t | n_dofs () const |
| std::size_t | 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 FEType *type=NULL, const std::set< std::string > *system_names=NULL) 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 std::set< std::string > *system_names=NULL) const |
| void | get_solution (std::vector< Number > &soln, std::vector< std::string > &names) 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=libMesh::out) |
| static unsigned int | n_objects () |
| static void | enable_print_counter_info () |
| static void | disable_print_counter_info () |
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 |
| static bool | _enable_print_counter = true |
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 67 of file equation_systems.h.
Member Typedef Documentation
typedef std::map<std::string, System*>::const_iterator libMesh::EquationSystems::const_system_iterator [protected] |
Typedef for constatnt system iterators
Definition at line 439 of file equation_systems.h.
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 113 of file reference_counter.h.
typedef std::map<std::string, System*>::iterator libMesh::EquationSystems::system_iterator [protected] |
Typedef for system iterators
Definition at line 434 of file equation_systems.h.
Member Enumeration Documentation
Define enumeration to set properties in EquationSystems::read()
- Enumerator:
READ_HEADER READ_DATA READ_ADDITIONAL_DATA READ_LEGACY_FORMAT TRY_READ_IFEMS READ_BASIC_ONLY
Definition at line 74 of file equation_systems.h.
00074 { READ_HEADER = 1, 00075 READ_DATA = 2, 00076 READ_ADDITIONAL_DATA = 4, 00077 READ_LEGACY_FORMAT = 8, 00078 TRY_READ_IFEMS = 16, 00079 READ_BASIC_ONLY = 32 };
Define enumeration to set properties in EquationSystems::write()
Definition at line 84 of file equation_systems.h.
00084 { WRITE_DATA = 1, 00085 WRITE_ADDITIONAL_DATA = 2, 00086 WRITE_PARALLEL_FILES = 4, 00087 WRITE_SERIAL_FILES = 8 };
Constructor & Destructor Documentation
Constructor.
Definition at line 54 of file equation_systems.C.
References parameters, libMesh::Real, libMesh::Parameters::set(), and libMesh::TOLERANCE.
00054 : 00055 _mesh (m), 00056 _mesh_data (mesh_data) 00057 { 00058 // Set default parameters 00059 this->parameters.set<Real> ("linear solver tolerance") = TOLERANCE * TOLERANCE; 00060 this->parameters.set<unsigned int>("linear solver maximum iterations") = 5000; 00061 }
| libMesh::EquationSystems::~EquationSystems | ( | ) | [virtual] |
Destructor. Should be virtual, since the user may want to derive subclasses of EquationSystems.
Definition at line 65 of file equation_systems.C.
References clear().
00066 { 00067 this->clear (); 00068 }
Member Function Documentation
| void libMesh::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 1232 of file equation_systems.C.
References _mesh, libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), libMesh::MeshBase::nodes_begin(), and libMesh::MeshBase::nodes_end().
Referenced by add_system().
01233 { 01234 // All the nodes 01235 MeshBase::node_iterator node_it = _mesh.nodes_begin(); 01236 const MeshBase::node_iterator node_end = _mesh.nodes_end(); 01237 01238 for ( ; node_it != node_end; ++node_it) 01239 (*node_it)->add_system(); 01240 01241 // All the elements 01242 MeshBase::element_iterator elem_it = _mesh.elements_begin(); 01243 const MeshBase::element_iterator elem_end = _mesh.elements_end(); 01244 01245 for ( ; elem_it != elem_end; ++elem_it) 01246 (*elem_it)->add_system(); 01247 }
| void libMesh::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 148 of file equation_systems_io.C.
References _mesh, add_system(), libMesh::Parallel::Communicator::broadcast(), libMesh::Xdr::close(), libMesh::CommWorld, libMesh::Xdr::data(), libMesh::MeshBase::fix_broken_node_and_element_numbering(), get_mesh(), get_system(), libMesh::MeshTools::Private::globally_renumber_nodes_and_elements(), init(), mesh, libMesh::processor_id(), read(), READ_ADDITIONAL_DATA, READ_BASIC_ONLY, READ_DATA, libMesh::System::read_header(), READ_HEADER, READ_LEGACY_FORMAT, libMesh::Xdr::reading(), libMesh::System::set_basic_system_only(), libMesh::Xdr::set_version(), TRY_READ_IFEMS, and update().
Referenced by read().
00151 { 00215 // Set booleans from the read_flags argument 00216 const bool read_header = read_flags & EquationSystems::READ_HEADER; 00217 const bool read_data = read_flags & EquationSystems::READ_DATA; 00218 const bool read_additional_data = read_flags & EquationSystems::READ_ADDITIONAL_DATA; 00219 const bool read_legacy_format = read_flags & EquationSystems::READ_LEGACY_FORMAT; 00220 const bool try_read_ifems = read_flags & EquationSystems::TRY_READ_IFEMS; 00221 const bool read_basic_only = read_flags & EquationSystems::READ_BASIC_ONLY; 00222 bool read_parallel_files = false; 00223 00224 std::map<std::string, System*> xda_systems; 00225 00226 // This will unzip a file with .bz2 as the extension, otherwise it 00227 // simply returns the name if the file need not be unzipped. 00228 Xdr io ((libMesh::processor_id() == 0) ? name : "", mode); 00229 libmesh_assert (io.reading()); 00230 00231 { 00232 // 1.) 00233 // Read the version header. 00234 std::string version = "legacy"; 00235 if (!read_legacy_format) 00236 { 00237 if (libMesh::processor_id() == 0) io.data(version); 00238 CommWorld.broadcast(version); 00239 00240 // All processors have the version header, if it does not contain 00241 // "libMesh" something then it is a legacy file. 00242 std::string::size_type lm_pos = version.find("libMesh"); 00243 if (!(lm_pos < version.size())) 00244 { 00245 io.close(); 00246 00247 // Recursively call this read() function but with the 00248 // EquationSystems::READ_LEGACY_FORMAT bit set. 00249 this->read (name, mode, (read_flags | EquationSystems::READ_LEGACY_FORMAT)); 00250 return; 00251 } 00252 00253 // Figure out the libMesh version that created this file 00254 std::istringstream iss(version.substr(lm_pos + 8)); 00255 int ver_major = 0, ver_minor = 0, ver_patch = 0; 00256 char dot; 00257 iss >> ver_major >> dot >> ver_minor >> dot >> ver_patch; 00258 io.set_version(LIBMESH_VERSION_ID(ver_major, ver_minor, ver_patch)); 00259 00260 00261 read_parallel_files = (version.rfind(" parallel") < version.size()); 00262 00263 // If requested that we try to read infinite element information, 00264 // and the string " with infinite elements" is not in the version, 00265 // then tack it on. This is for compatibility reading ifem 00266 // files written prior to 11/10/2008 - BSK 00267 if (try_read_ifems) 00268 if (!(version.rfind(" with infinite elements") < version.size())) 00269 version += " with infinite elements"; 00270 00271 } 00272 else 00273 libmesh_deprecated(); 00274 00275 START_LOG("read()","EquationSystems"); 00276 00277 // 2.) 00278 // Read the number of equation systems 00279 unsigned int n_sys=0; 00280 if (libMesh::processor_id() == 0) io.data (n_sys); 00281 CommWorld.broadcast(n_sys); 00282 00283 for (unsigned int sys=0; sys<n_sys; sys++) 00284 { 00285 // 3.) 00286 // Read the name of the sys-th equation system 00287 std::string sys_name; 00288 if (libMesh::processor_id() == 0) io.data (sys_name); 00289 CommWorld.broadcast(sys_name); 00290 00291 // 4.) 00292 // Read the type of the sys-th equation system 00293 std::string sys_type; 00294 if (libMesh::processor_id() == 0) io.data (sys_type); 00295 CommWorld.broadcast(sys_type); 00296 00297 if (read_header) 00298 this->add_system (sys_type, sys_name); 00299 00300 // 5.) - 9.) 00301 // Let System::read_header() do the job 00302 System& new_system = this->get_system(sys_name); 00303 new_system.read_header (io, 00304 version, 00305 read_header, 00306 read_additional_data, 00307 read_legacy_format); 00308 00309 xda_systems.insert(std::make_pair(sys_name, &new_system)); 00310 00311 // If we're only creating "basic" systems, we need to tell 00312 // each system that before we call init() later. 00313 if (read_basic_only) 00314 new_system.set_basic_system_only(); 00315 } 00316 } 00317 00318 00319 00320 // Now we are ready to initialize the underlying data 00321 // structures. This will initialize the vectors for 00322 // storage, the dof_map, etc... 00323 if (read_header) 00324 this->init(); 00325 00326 // 10.) & 11.) 00327 // Read and set the numeric vector values 00328 if (read_data) 00329 { 00330 // the EquationSystems::read() method should look constant from the mesh 00331 // perspective, but we need to assign a temporary numbering to the nodes 00332 // and elements in the mesh, which requires that we abuse const_cast 00333 if (!read_legacy_format) 00334 { 00335 MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh()); 00336 MeshTools::Private::globally_renumber_nodes_and_elements(mesh); 00337 } 00338 00339 Xdr local_io (read_parallel_files ? local_file_name(name) : "", mode); 00340 00341 std::map<std::string, System*>::iterator 00342 pos = xda_systems.begin(); 00343 00344 for (; pos != xda_systems.end(); ++pos) 00345 if (read_legacy_format) 00346 { 00347 libmesh_deprecated(); 00348 pos->second->read_legacy_data (io, read_additional_data); 00349 } 00350 else 00351 if (read_parallel_files) 00352 pos->second->read_parallel_data (local_io, read_additional_data); 00353 else 00354 pos->second->read_serialized_data (io, read_additional_data); 00355 00356 00357 // Undo the temporary numbering. 00358 if (!read_legacy_format) 00359 _mesh.fix_broken_node_and_element_numbering(); 00360 } 00361 00362 STOP_LOG("read()","EquationSystems"); 00363 00364 // Localize each system's data 00365 this->update(); 00366 }
| T_sys & libMesh::EquationSystems::add_system | ( | const std::string & | name | ) | [inline] |
Add the system named name to the systems array.
Definition at line 512 of file equation_systems.h.
References _add_system_to_nodes_and_elems(), _systems, and n_systems().
00513 { 00514 T_sys* ptr = NULL; 00515 00516 if (!_systems.count(name)) 00517 { 00518 ptr = new T_sys(*this, name, this->n_systems()); 00519 00520 _systems.insert (std::make_pair(name, ptr)); 00521 00522 // Tell all the \p DofObject entities to add a system. 00523 this->_add_system_to_nodes_and_elems(); 00524 } 00525 else 00526 { 00527 // We now allow redundant add_system calls, to make it 00528 // easier to load data from files for user-derived system 00529 // subclasses 00530 // libMesh::err << "ERROR: There was already a system" 00531 // << " named " << name 00532 // << std::endl; 00533 00534 // libmesh_error(); 00535 00536 ptr = &(this->get_system<T_sys>(name)); 00537 } 00538 00539 // Return a dynamically casted reference to the newly added System. 00540 return *ptr; 00541 }
| System & libMesh::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 335 of file equation_systems.C.
References _systems, libMesh::err, and get_system().
Referenced by _read_impl(), and libMesh::ErrorVector::plot_error().
00337 { 00338 // If the user already built a system with this name, we'll 00339 // trust them and we'll use it. That way they can pre-add 00340 // non-standard derived system classes, and if their restart file 00341 // has some non-standard sys_type we won't throw an error. 00342 if (_systems.count(name)) 00343 { 00344 return this->get_system(name); 00345 } 00346 // Build a basic System 00347 else if (sys_type == "Basic") 00348 this->add_system<System> (name); 00349 00350 // Build a Newmark system 00351 else if (sys_type == "Newmark") 00352 this->add_system<NewmarkSystem> (name); 00353 00354 // Build an Explicit system 00355 else if ((sys_type == "Explicit")) 00356 this->add_system<ExplicitSystem> (name); 00357 00358 // Build an Implicit system 00359 else if ((sys_type == "Implicit") || 00360 (sys_type == "Steady" )) 00361 this->add_system<ImplicitSystem> (name); 00362 00363 // build a transient implicit linear system 00364 else if ((sys_type == "Transient") || 00365 (sys_type == "TransientImplicit") || 00366 (sys_type == "TransientLinearImplicit")) 00367 this->add_system<TransientLinearImplicitSystem> (name); 00368 00369 // build a transient implicit nonlinear system 00370 else if (sys_type == "TransientNonlinearImplicit") 00371 this->add_system<TransientNonlinearImplicitSystem> (name); 00372 00373 // build a transient explicit system 00374 else if (sys_type == "TransientExplicit") 00375 this->add_system<TransientExplicitSystem> (name); 00376 00377 // build a linear implicit system 00378 else if (sys_type == "LinearImplicit") 00379 this->add_system<LinearImplicitSystem> (name); 00380 00381 // build a nonlinear implicit system 00382 else if (sys_type == "NonlinearImplicit") 00383 this->add_system<NonlinearImplicitSystem> (name); 00384 00385 // build a Reduced Basis Construction system 00386 else if (sys_type == "RBConstruction") 00387 this->add_system<RBConstruction> (name); 00388 00389 // build a transient Reduced Basis Construction system 00390 else if (sys_type == "TransientRBConstruction") 00391 this->add_system<TransientRBConstruction> (name); 00392 00393 #ifdef LIBMESH_HAVE_SLEPC 00394 // build an eigen system 00395 else if (sys_type == "Eigen") 00396 this->add_system<EigenSystem> (name); 00397 #endif 00398 00399 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) 00400 // build a frequency system 00401 else if (sys_type == "Frequency") 00402 this->add_system<FrequencySystem> (name); 00403 #endif 00404 00405 else 00406 { 00407 libMesh::err << "ERROR: Unknown system type: " 00408 << sys_type 00409 << std::endl; 00410 libmesh_error(); 00411 } 00412 00413 // Return a reference to the new system 00414 //return (*this)(name); 00415 return this->get_system(name); 00416 }
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 462 of file equation_systems.C.
References get_system(), and n_systems().
00463 { 00464 libmesh_assert (this->n_systems()); 00465 00466 for (unsigned int i=this->n_systems(); i != 0; --i) 00467 this->get_system(i-1).adjoint_solve(qoi_indices); 00468 }
| void libMesh::EquationSystems::allgather | ( | ) |
Serializes a distributed mesh and its associated degree of freedom numbering for all systems
Definition at line 271 of file equation_systems.C.
References _mesh, libMesh::MeshBase::allgather(), libMesh::DofMap::create_dof_constraints(), libMesh::DofMap::distribute_dofs(), libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), libMesh::System::get_dof_map(), get_system(), libMesh::MeshBase::is_serial(), n_systems(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), libMesh::DofMap::prepare_send_list(), libMesh::DofMap::process_constraints(), libMesh::System::time, and libMesh::System::user_constrain().
00272 { 00273 // A serial mesh means nothing needs to be done 00274 if (_mesh.is_serial()) 00275 return; 00276 00277 const unsigned int n_sys = this->n_systems(); 00278 00279 libmesh_assert_not_equal_to (n_sys, 0); 00280 00281 // Gather the mesh 00282 _mesh.allgather(); 00283 00284 // Tell all the \p DofObject entities how many systems 00285 // there are. 00286 { 00287 MeshBase::node_iterator node_it = _mesh.nodes_begin(); 00288 const MeshBase::node_iterator node_end = _mesh.nodes_end(); 00289 00290 for ( ; node_it != node_end; ++node_it) 00291 (*node_it)->set_n_systems(n_sys); 00292 00293 MeshBase::element_iterator elem_it = _mesh.elements_begin(); 00294 const MeshBase::element_iterator elem_end = _mesh.elements_end(); 00295 00296 for ( ; elem_it != elem_end; ++elem_it) 00297 (*elem_it)->set_n_systems(n_sys); 00298 } 00299 00300 // And distribute each system's dofs 00301 for (unsigned int i=0; i != this->n_systems(); ++i) 00302 { 00303 System &sys = this->get_system(i); 00304 DofMap &dof_map = sys.get_dof_map(); 00305 dof_map.distribute_dofs(_mesh); 00306 00307 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00308 // The user probably won't need constraint equations or the 00309 // send_list after an allgather, but let's keep it in consistent 00310 // shape just in case. 00311 dof_map.create_dof_constraints(_mesh, sys.time); 00312 sys.user_constrain(); 00313 dof_map.process_constraints(_mesh); 00314 #endif 00315 dof_map.prepare_send_list(); 00316 } 00317 }
| void libMesh::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 971 of file equation_systems.C.
References _mesh, _systems, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::DofMap::dof_indices(), end, libMesh::System::get_dof_map(), libMesh::Elem::infinite(), libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_nodes(), n_systems(), libMesh::System::n_vars(), n_vars(), libMesh::FEInterface::nodal_soln(), libMesh::MeshBase::processor_id(), libMesh::System::update_global_solution(), and libMesh::System::variable_type().
Referenced by libMesh::ExodusII_IO::write_discontinuous_exodusII(), and libMesh::GMVIO::write_discontinuous_gmv().
00972 { 00973 START_LOG("build_discontinuous_solution_vector()", "EquationSystems"); 00974 00975 libmesh_assert (this->n_systems()); 00976 00977 const unsigned int dim = _mesh.mesh_dimension(); 00978 const unsigned int nv = this->n_vars(); 00979 unsigned int tw=0; 00980 00981 // get the total weight 00982 { 00983 MeshBase::element_iterator it = _mesh.active_elements_begin(); 00984 const MeshBase::element_iterator end = _mesh.active_elements_end(); 00985 00986 for ( ; it != end; ++it) 00987 tw += (*it)->n_nodes(); 00988 } 00989 00990 00991 // Only if we are on processor zero, allocate the storage 00992 // to hold (number_of_nodes)*(number_of_variables) entries. 00993 if (_mesh.processor_id() == 0) 00994 soln.resize(tw*nv); 00995 00996 std::vector<Number> sys_soln; 00997 00998 00999 unsigned int var_num=0; 01000 01001 // For each system in this EquationSystems object, 01002 // update the global solution and if we are on processor 0, 01003 // loop over the elements and build the nodal solution 01004 // from the element solution. Then insert this nodal solution 01005 // into the vector passed to build_solution_vector. 01006 const_system_iterator pos = _systems.begin(); 01007 const const_system_iterator end = _systems.end(); 01008 01009 for (; pos != end; ++pos) 01010 { 01011 const System& system = *(pos->second); 01012 const unsigned int nv_sys = system.n_vars(); 01013 01014 system.update_global_solution (sys_soln, 0); 01015 01016 if (_mesh.processor_id() == 0) 01017 { 01018 std::vector<Number> elem_soln; // The finite element solution 01019 std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes 01020 std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element 01021 01022 for (unsigned int var=0; var<nv_sys; var++) 01023 { 01024 const FEType& fe_type = system.variable_type(var); 01025 01026 MeshBase::element_iterator it = _mesh.active_elements_begin(); 01027 const MeshBase::element_iterator end_elem = _mesh.active_elements_end(); 01028 01029 unsigned int nn=0; 01030 01031 for ( ; it != end_elem; ++it) 01032 { 01033 const Elem* elem = *it; 01034 system.get_dof_map().dof_indices (elem, dof_indices, var); 01035 01036 elem_soln.resize(dof_indices.size()); 01037 01038 for (unsigned int i=0; i<dof_indices.size(); i++) 01039 elem_soln[i] = sys_soln[dof_indices[i]]; 01040 01041 FEInterface::nodal_soln (dim, 01042 fe_type, 01043 elem, 01044 elem_soln, 01045 nodal_soln); 01046 01047 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 01048 // infinite elements should be skipped... 01049 if (!elem->infinite()) 01050 #endif 01051 { 01052 libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes()); 01053 01054 for (unsigned int n=0; n<elem->n_nodes(); n++) 01055 { 01056 soln[nv*(nn++) + (var + var_num)] += 01057 nodal_soln[n]; 01058 } 01059 } 01060 } 01061 } 01062 } 01063 01064 var_num += nv_sys; 01065 } 01066 01067 STOP_LOG("build_discontinuous_solution_vector()", "EquationSystems"); 01068 }
| void libMesh::EquationSystems::build_solution_vector | ( | std::vector< Number > & | soln, | |
| const std::set< std::string > * | system_names = NULL | |||
| ) | 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()). If systems_names!=NULL, only include data from the specified systems.
Definition at line 673 of file equation_systems.C.
References _mesh, _systems, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::Variable::active_on_subdomain(), libMesh::CommWorld, libMesh::DofMap::dof_indices(), end, libMesh::FEInterface::field_type(), libMesh::System::get_dof_map(), libMesh::Variable::implicitly_active(), libMesh::Elem::infinite(), std::max(), libMesh::MeshBase::max_node_id(), libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_nodes(), libMesh::MeshBase::n_nodes(), n_systems(), libMesh::System::n_vars(), libMesh::FEInterface::n_vec_dim(), libMesh::FEInterface::nodal_soln(), libMesh::Elem::node(), libMesh::Real, libMesh::Parallel::Communicator::sum(), libMeshEnums::TYPE_VECTOR, libMesh::System::update_global_solution(), libMesh::System::variable(), libMesh::System::variable_type(), and libMesh::zero.
00675 { 00676 START_LOG("build_solution_vector()", "EquationSystems"); 00677 00678 // This function must be run on all processors at once 00679 parallel_only(); 00680 00681 libmesh_assert (this->n_systems()); 00682 00683 const unsigned int dim = _mesh.mesh_dimension(); 00684 const dof_id_type nn = _mesh.n_nodes(); 00685 const unsigned short int one = 1; 00686 00687 // We'd better have a contiguous node numbering 00688 libmesh_assert_equal_to (nn, _mesh.max_node_id()); 00689 00690 // allocate storage to hold 00691 // (number_of_nodes)*(number_of_variables) entries. 00692 // We have to differentiate between between scalar and vector 00693 // variables. We intercept vector variables and treat each 00694 // component as a scalar variable (consistently with build_solution_names). 00695 00696 unsigned int nv = 0; 00697 00698 //Could this be replaced by a/some convenience methods?[PB] 00699 { 00700 unsigned int n_scalar_vars = 0; 00701 unsigned int n_vector_vars = 0; 00702 const_system_iterator pos = _systems.begin(); 00703 const const_system_iterator end = _systems.end(); 00704 00705 for (; pos != end; ++pos) 00706 { 00707 // Check current system is listed in system_names, and skip pos if not 00708 bool use_current_system = (system_names == NULL); 00709 if(!use_current_system) 00710 { 00711 use_current_system = std::find( system_names->begin(), system_names->end(), pos->first ) != system_names->end(); 00712 } 00713 if(!use_current_system) 00714 { 00715 continue; 00716 } 00717 00718 for (unsigned int vn=0; vn<pos->second->n_vars(); vn++) 00719 { 00720 if( FEInterface::field_type(pos->second->variable_type(vn)) == 00721 TYPE_VECTOR ) 00722 n_vector_vars++; 00723 else 00724 n_scalar_vars++; 00725 } 00726 } 00727 // Here, we're assuming the number of vector components is the same 00728 // as the mesh dimension. Will break for mixed dimension meshes. 00729 nv = n_scalar_vars + dim*n_vector_vars; 00730 } 00731 00732 soln.resize(nn*nv); 00733 00734 // Zero out the soln vector 00735 std::fill (soln.begin(), soln.end(), libMesh::zero); 00736 00737 std::vector<Number> sys_soln; 00738 00739 // (Note that we use an unsigned short int here even though an 00740 // unsigned char would be more that sufficient. The MPI 1.1 00741 // standard does not require that MPI_SUM, MPI_PROD etc... be 00742 // implemented for char data types. 12/23/2003 - BSK) 00743 std::vector<unsigned short int> node_conn(nn), repeat_count(nn); 00744 00745 // Get the number of elements that share each node. We will 00746 // compute the average value at each node. This is particularly 00747 // useful for plotting discontinuous data. 00748 MeshBase::element_iterator e_it = _mesh.active_local_elements_begin(); 00749 const MeshBase::element_iterator e_end = _mesh.active_local_elements_end(); 00750 00751 for ( ; e_it != e_end; ++e_it) 00752 for (unsigned int n=0; n<(*e_it)->n_nodes(); n++) 00753 node_conn[(*e_it)->node(n)]++; 00754 00755 // Gather the distributed node_conn arrays in the case of 00756 // multiple processors. 00757 CommWorld.sum(node_conn); 00758 00759 unsigned int var_num=0; 00760 00761 // For each system in this EquationSystems object, 00762 // update the global solution and if we are on processor 0, 00763 // loop over the elements and build the nodal solution 00764 // from the element solution. Then insert this nodal solution 00765 // into the vector passed to build_solution_vector. 00766 const_system_iterator pos = _systems.begin(); 00767 const const_system_iterator end = _systems.end(); 00768 00769 for (; pos != end; ++pos) 00770 { 00771 // Check current system is listed in system_names, and skip pos if not 00772 bool use_current_system = (system_names == NULL); 00773 if(!use_current_system) 00774 { 00775 use_current_system = std::find( system_names->begin(), system_names->end(), pos->first ) != system_names->end(); 00776 } 00777 if(!use_current_system) 00778 { 00779 continue; 00780 } 00781 00782 const System& system = *(pos->second); 00783 const unsigned int nv_sys = system.n_vars(); 00784 00785 //Could this be replaced by a/some convenience methods?[PB] 00786 unsigned int n_scalar_vars = 0; 00787 unsigned int n_vector_vars = 0; 00788 for (unsigned int vn=0; vn<pos->second->n_vars(); vn++) 00789 { 00790 if( FEInterface::field_type(pos->second->variable_type(vn)) == 00791 TYPE_VECTOR ) 00792 n_vector_vars++; 00793 else 00794 n_scalar_vars++; 00795 } 00796 00797 // Here, we're assuming the number of vector components is the same 00798 // as the mesh dimension. Will break for mixed dimension meshes. 00799 unsigned int nv_sys_split = n_scalar_vars + dim*n_vector_vars; 00800 00801 system.update_global_solution (sys_soln); 00802 00803 std::vector<Number> elem_soln; // The finite element solution 00804 std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes 00805 std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element 00806 00807 for (unsigned int var=0; var<nv_sys; var++) 00808 { 00809 const FEType& fe_type = system.variable_type(var); 00810 const Variable &var_description = system.variable(var); 00811 const DofMap &dof_map = system.get_dof_map(); 00812 00813 unsigned int n_vec_dim = FEInterface::n_vec_dim( pos->second->get_mesh(), fe_type ); 00814 00815 std::fill (repeat_count.begin(), repeat_count.end(), 0); 00816 00817 MeshBase::element_iterator it = _mesh.active_local_elements_begin(); 00818 const MeshBase::element_iterator end_elem = _mesh.active_local_elements_end(); 00819 00820 for ( ; it != end_elem; ++it) 00821 if (var_description.active_on_subdomain((*it)->subdomain_id())) 00822 { 00823 const Elem* elem = *it; 00824 00825 dof_map.dof_indices (elem, dof_indices, var); 00826 00827 elem_soln.resize(dof_indices.size()); 00828 00829 for (unsigned int i=0; i<dof_indices.size(); i++) 00830 elem_soln[i] = sys_soln[dof_indices[i]]; 00831 00832 FEInterface::nodal_soln (dim, 00833 fe_type, 00834 elem, 00835 elem_soln, 00836 nodal_soln); 00837 00838 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00839 // infinite elements should be skipped... 00840 if (!elem->infinite()) 00841 #endif 00842 { 00843 libmesh_assert_equal_to (nodal_soln.size(), n_vec_dim*elem->n_nodes()); 00844 00845 for (unsigned int n=0; n<elem->n_nodes(); n++) 00846 { 00847 repeat_count[elem->node(n)]++; 00848 for( unsigned int d=0; d < n_vec_dim; d++ ) 00849 { 00850 // For vector-valued elements, all components are in nodal_soln. For each 00851 // node, the components are stored in order, i.e. node_0 -> s0_x, s0_y, s0_z 00852 soln[nv*(elem->node(n)) + (var+d + var_num)] += nodal_soln[n_vec_dim*n+d]; 00853 } 00854 } 00855 } 00856 } // end loop over elements 00857 00858 // when a variable is active everywhere the repeat_count 00859 // and node_conn arrays should be identical, so let's 00860 // use that information to avoid unnecessary communication 00861 if (var_description.implicitly_active()) 00862 repeat_count = node_conn; 00863 00864 else 00865 CommWorld.sum (repeat_count); 00866 00867 for (unsigned int n=0; n<nn; n++) 00868 for( unsigned int d=0; d < n_vec_dim; d++ ) 00869 soln[nv*n + (var+d + var_num)] /= 00870 static_cast<Real>(std::max (repeat_count[n], one)); 00871 00872 } // end loop on variables in this system 00873 00874 var_num += nv_sys_split; 00875 } // end loop over systems 00876 00877 // Now each processor has computed contriburions to the 00878 // soln vector. Gather them all up. 00879 CommWorld.sum(soln); 00880 00881 STOP_LOG("build_solution_vector()", "EquationSystems"); 00882 }
| void libMesh::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 586 of file equation_systems.C.
Referenced by libMesh::MeshOutput< MT >::_build_variable_names_and_solution_vector().
00589 { 00590 //TODO:[BSK] re-implement this from the method below 00591 libmesh_error(); 00592 00593 // // Get a reference to the named system 00594 // const System& system = this->get_system(system_name); 00595 00596 // // Get the number associated with the variable_name we are passed 00597 // const unsigned short int variable_num = system.variable_number(variable_name); 00598 00599 // // Get the dimension of the current mesh 00600 // const unsigned int dim = _mesh.mesh_dimension(); 00601 00602 // // If we're on processor 0, allocate enough memory to hold the solution. 00603 // // Since we're only looking at one variable, there will be one solution value 00604 // // for each node in the mesh. 00605 // if (_mesh.processor_id() == 0) 00606 // soln.resize(_mesh.n_nodes()); 00607 00608 // // Vector to hold the global solution from all processors 00609 // std::vector<Number> sys_soln; 00610 00611 // // Update the global solution from all processors 00612 // system.update_global_solution (sys_soln, 0); 00613 00614 // // Temporary vector to store the solution on an individual element. 00615 // std::vector<Number> elem_soln; 00616 00617 // // The FE solution interpolated to the nodes 00618 // std::vector<Number> nodal_soln; 00619 00620 // // The DOF indices for the element 00621 // std::vector<dof_id_type> dof_indices; 00622 00623 // // Determine the finite/infinite element type used in this system 00624 // const FEType& fe_type = system.variable_type(variable_num); 00625 00626 // // Define iterators to iterate over all the elements of the mesh 00627 // const_active_elem_iterator it (_mesh.elements_begin()); 00628 // const const_active_elem_iterator end(_mesh.elements_end()); 00629 00630 // // Loop over elements 00631 // for ( ; it != end; ++it) 00632 // { 00633 // // Convenient shortcut to the element pointer 00634 // const Elem* elem = *it; 00635 00636 // // Fill the dof_indices vector for this variable 00637 // system.get_dof_map().dof_indices(elem, 00638 // dof_indices, 00639 // variable_num); 00640 00641 // // Resize the element solution vector to fit the 00642 // // dof_indices for this element. 00643 // elem_soln.resize(dof_indices.size()); 00644 00645 // // Transfer the system solution to the element 00646 // // solution by mapping it through the dof_indices vector. 00647 // for (unsigned int i=0; i<dof_indices.size(); i++) 00648 // elem_soln[i] = sys_soln[dof_indices[i]]; 00649 00650 // // Using the FE interface, compute the nodal_soln 00651 // // for the current elemnt type given the elem_soln 00652 // FEInterface::nodal_soln (dim, 00653 // fe_type, 00654 // elem, 00655 // elem_soln, 00656 // nodal_soln); 00657 00658 // // Sanity check -- make sure that there are the same number 00659 // // of entries in the nodal_soln as there are nodes in the 00660 // // element! 00661 // libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes()); 00662 00663 // // Copy the nodal solution over into the correct place in 00664 // // the global soln vector which will be returned to the user. 00665 // for (unsigned int n=0; n<elem->n_nodes(); n++) 00666 // soln[elem->node(n)] = nodal_soln[n]; 00667 // } 00668 }
| void libMesh::EquationSystems::build_variable_names | ( | std::vector< std::string > & | var_names, | |
| const FEType * | type = NULL, |
|||
| const std::set< std::string > * | system_names = NULL | |||
| ) | const |
Fill the input vector var_names with the names of the variables for each system. If type is passed, only variables of the specified type will be populated. If systems_names!=NULL, only include names from the specified systems.
Definition at line 472 of file equation_systems.C.
References _systems, end, libMesh::FEInterface::field_type(), get_mesh(), libMesh::MeshBase::mesh_dimension(), n_systems(), n_vars(), libMesh::FEInterface::n_vec_dim(), and libMeshEnums::TYPE_VECTOR.
Referenced by libMesh::MeshOutput< MT >::_build_variable_names_and_solution_vector(), libMesh::ExodusII_IO::write_discontinuous_exodusII(), and libMesh::GMVIO::write_discontinuous_gmv().
00475 { 00476 libmesh_assert (this->n_systems()); 00477 00478 unsigned int var_num=0; 00479 00480 const_system_iterator pos = _systems.begin(); 00481 const const_system_iterator end = _systems.end(); 00482 00483 // Need to size var_names by scalar variables plus all the 00484 // vector components for all the vector variables 00485 //Could this be replaced by a/some convenience methods?[PB] 00486 { 00487 unsigned int n_scalar_vars = 0; 00488 unsigned int n_vector_vars = 0; 00489 00490 for (; pos != end; ++pos) 00491 { 00492 // Check current system is listed in system_names, and skip pos if not 00493 bool use_current_system = (system_names == NULL); 00494 if(!use_current_system) 00495 { 00496 use_current_system = std::find( system_names->begin(), system_names->end(), pos->first ) != system_names->end(); 00497 } 00498 if(!use_current_system) 00499 { 00500 continue; 00501 } 00502 00503 for (unsigned int vn=0; vn<pos->second->n_vars(); vn++) 00504 { 00505 if( FEInterface::field_type(pos->second->variable_type(vn)) == 00506 TYPE_VECTOR ) 00507 n_vector_vars++; 00508 else 00509 n_scalar_vars++; 00510 } 00511 } 00512 00513 // Here, we're assuming the number of vector components is the same 00514 // as the mesh dimension. Will break for mixed dimension meshes. 00515 unsigned int dim = this->get_mesh().mesh_dimension(); 00516 unsigned int nv = n_scalar_vars + dim*n_vector_vars; 00517 00518 // We'd better not have more than dim*his->n_vars() (all vector variables) 00519 libmesh_assert_less_equal ( nv, dim*this->n_vars() ); 00520 00521 // Here, we're assuming the number of vector components is the same 00522 // as the mesh dimension. Will break for mixed dimension meshes. 00523 00524 var_names.resize( nv ); 00525 } 00526 00527 // reset 00528 pos = _systems.begin(); 00529 00530 for (; pos != end; ++pos) 00531 { 00532 // Check current system is listed in system_names, and skip pos if not 00533 bool use_current_system = (system_names == NULL); 00534 if(!use_current_system) 00535 { 00536 use_current_system = std::find( system_names->begin(), system_names->end(), pos->first ) != system_names->end(); 00537 } 00538 if(!use_current_system) 00539 { 00540 continue; 00541 } 00542 00543 for (unsigned int vn=0; vn<pos->second->n_vars(); vn++) 00544 { 00545 std::string var_name = pos->second->variable_name(vn); 00546 FEType fe_type = pos->second->variable_type(vn); 00547 00548 unsigned int n_vec_dim = FEInterface::n_vec_dim( pos->second->get_mesh(), fe_type); 00549 00550 // Filter on the type if requested 00551 if (type == NULL || (type && *type == fe_type)) 00552 { 00553 if( FEInterface::field_type(fe_type) == TYPE_VECTOR ) 00554 { 00555 switch(n_vec_dim) 00556 { 00557 case 0: 00558 case 1: 00559 var_names[var_num++] = var_name; 00560 break; 00561 case 2: 00562 var_names[var_num++] = var_name+"_x"; 00563 var_names[var_num++] = var_name+"_y"; 00564 break; 00565 case 3: 00566 var_names[var_num++] = var_name+"_x"; 00567 var_names[var_num++] = var_name+"_y"; 00568 var_names[var_num++] = var_name+"_z"; 00569 break; 00570 default: 00571 std::cerr << "Invalid dim in build_variable_names" << std::endl; 00572 libmesh_error(); 00573 } 00574 } 00575 else 00576 var_names[var_num++] = var_name; 00577 } 00578 } 00579 } 00580 // Now resize again in case we filtered any names 00581 var_names.resize(var_num); 00582 }
| void libMesh::EquationSystems::clear | ( | ) | [virtual] |
Returns tha data structure to a pristine state.
Definition at line 72 of file equation_systems.C.
References _systems, libMesh::Parameters::clear(), and parameters.
Referenced by ~EquationSystems().
00073 { 00074 // Clear any additional parameters 00075 parameters.clear (); 00076 00077 // clear the systems. We must delete them 00078 // since we newed them! 00079 while (!_systems.empty()) 00080 { 00081 system_iterator pos = _systems.begin(); 00082 00083 System *sys = pos->second; 00084 delete sys; 00085 sys = NULL; 00086 00087 _systems.erase (pos); 00088 } 00089 }
| bool libMesh::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 1072 of file equation_systems.C.
References _systems, libMesh::System::compare(), end, get_system(), n_systems(), and libMesh::out.
01075 { 01076 // safety check, whether we handle at least the same number 01077 // of systems 01078 std::vector<bool> os_result; 01079 01080 if (this->n_systems() != other_es.n_systems()) 01081 { 01082 if (verbose) 01083 { 01084 libMesh::out << " Fatal difference. This system handles " 01085 << this->n_systems() << " systems," << std::endl 01086 << " while the other system handles " 01087 << other_es.n_systems() 01088 << " systems." << std::endl 01089 << " Aborting comparison." << std::endl; 01090 } 01091 return false; 01092 } 01093 else 01094 { 01095 // start comparing each system 01096 const_system_iterator pos = _systems.begin(); 01097 const const_system_iterator end = _systems.end(); 01098 01099 for (; pos != end; ++pos) 01100 { 01101 const std::string& sys_name = pos->first; 01102 const System& system = *(pos->second); 01103 01104 // get the other system 01105 const System& other_system = other_es.get_system (sys_name); 01106 01107 os_result.push_back (system.compare (other_system, threshold, verbose)); 01108 01109 } 01110 01111 } 01112 01113 01114 // sum up the results 01115 if (os_result.size()==0) 01116 return true; 01117 else 01118 { 01119 bool os_identical; 01120 unsigned int n = 0; 01121 do 01122 { 01123 os_identical = os_result[n]; 01124 n++; 01125 } 01126 while (os_identical && n<os_result.size()); 01127 return os_identical; 01128 } 01129 }
| void libMesh::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 423 of file equation_systems.C.
References _systems, and libMesh::err.
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00107 { 00108 _enable_print_counter = false; 00109 return; 00110 }
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00101 { 00102 _enable_print_counter = true; 00103 return; 00104 }
| std::string libMesh::ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
00048 { 00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00050 00051 std::ostringstream oss; 00052 00053 oss << '\n' 00054 << " ---------------------------------------------------------------------------- \n" 00055 << "| Reference count information |\n" 00056 << " ---------------------------------------------------------------------------- \n"; 00057 00058 for (Counts::iterator it = _counts.begin(); 00059 it != _counts.end(); ++it) 00060 { 00061 const std::string name(it->first); 00062 const unsigned int creations = it->second.first; 00063 const unsigned int destructions = it->second.second; 00064 00065 oss << "| " << name << " reference count information:\n" 00066 << "| Creations: " << creations << '\n' 00067 << "| Destructions: " << destructions << '\n'; 00068 } 00069 00070 oss << " ---------------------------------------------------------------------------- \n"; 00071 00072 return oss.str(); 00073 00074 #else 00075 00076 return ""; 00077 00078 #endif 00079 }
| std::string libMesh::EquationSystems::get_info | ( | ) | const [virtual] |
- Returns:
- a string containing information about the systems, flags, and parameters.
Definition at line 1133 of file equation_systems.C.
References _systems, end, and n_systems().
Referenced by print_info().
01134 { 01135 std::ostringstream oss; 01136 01137 oss << " EquationSystems\n" 01138 << " n_systems()=" << this->n_systems() << '\n'; 01139 01140 // Print the info for the individual systems 01141 const_system_iterator pos = _systems.begin(); 01142 const const_system_iterator end = _systems.end(); 01143 01144 for (; pos != end; ++pos) 01145 oss << pos->second->get_info(); 01146 01147 01148 // // Possibly print the parameters 01149 // if (!this->parameters.empty()) 01150 // { 01151 // oss << " n_parameters()=" << this->n_parameters() << '\n'; 01152 // oss << " Parameters:\n"; 01153 01154 // for (std::map<std::string, Real>::const_iterator 01155 // param = _parameters.begin(); param != _parameters.end(); 01156 // ++param) 01157 // oss << " " 01158 // << "\"" 01159 // << param->first 01160 // << "\"" 01161 // << "=" 01162 // << param->second 01163 // << '\n'; 01164 // } 01165 01166 return oss.str(); 01167 }
| MeshBase & libMesh::EquationSystems::get_mesh | ( | ) | [inline] |
- Returns:
- a reference to the mesh
Definition at line 473 of file equation_systems.h.
References _mesh.
00474 { 00475 return _mesh; 00476 }
| const MeshBase & libMesh::EquationSystems::get_mesh | ( | ) | const [inline] |
- Returns:
- a constant reference to the mesh
Definition at line 465 of file equation_systems.h.
References _mesh.
Referenced by _read_impl(), build_variable_names(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), libMesh::MeshFunction::init(), libMesh::MeshFunction::operator()(), reinit(), libMesh::VTKIO::system_vectors_to_vtk(), and write().
00466 { 00467 return _mesh; 00468 }
| MeshData & libMesh::EquationSystems::get_mesh_data | ( | ) | [inline] |
- Returns:
- a reference to the mesh_data
Definition at line 488 of file equation_systems.h.
References _mesh_data.
00489 { 00490 libmesh_assert(_mesh_data); 00491 return *_mesh_data; 00492 }
| const MeshData & libMesh::EquationSystems::get_mesh_data | ( | ) | const [inline] |
- Returns:
- a constant reference to the mesh_data
Definition at line 480 of file equation_systems.h.
References _mesh_data.
00481 { 00482 libmesh_assert(_mesh_data); 00483 return *_mesh_data; 00484 }
| void libMesh::EquationSystems::get_solution | ( | std::vector< Number > & | soln, | |
| std::vector< std::string > & | names | |||
| ) | const |
Retrieve the solution data for CONSTANT MONOMIALs. If names is populated, only the variables corresponding to those names will be retrieved. This can be used to filter which variables are retrieved.
Definition at line 885 of file equation_systems.C.
References _mesh, _systems, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::Variable::active_on_subdomain(), libMesh::CommWorld, libMeshEnums::CONSTANT, libMesh::DofMap::dof_indices(), end, libMesh::System::get_dof_map(), libMesh::DofObject::id(), libMesh::MeshBase::max_elem_id(), libMeshEnums::MONOMIAL, libMesh::MeshBase::n_elem(), n_systems(), libMesh::System::n_vars(), libMesh::Parallel::Communicator::sum(), libMesh::System::update_global_solution(), libMesh::System::variable(), libMesh::System::variable_name(), and libMesh::System::variable_type().
00887 { 00888 // This function must be run on all processors at once 00889 parallel_only(); 00890 00891 libmesh_assert (this->n_systems()); 00892 00893 const dof_id_type ne = _mesh.n_elem(); 00894 00895 libmesh_assert_equal_to (ne, _mesh.max_elem_id()); 00896 00897 // If the names vector has entries, we will only populate the soln vector 00898 // with names included in that list. Note: The names vector may be 00899 // reordered upon exiting this function 00900 std::vector<std::string> filter_names = names; 00901 bool is_filter_names = ! filter_names.empty(); 00902 00903 soln.clear(); 00904 names.clear(); 00905 00906 std::vector<Number> sys_soln; 00907 00908 const FEType type(CONSTANT, MONOMIAL); 00909 00910 // For each system in this EquationSystems object, 00911 // update the global solution and collect the 00912 // CONSTANT MONOMIALs. The entries are in variable-major 00913 // format. 00914 const_system_iterator pos = _systems.begin(); 00915 const const_system_iterator end = _systems.end(); 00916 00917 for (; pos != end; ++pos) 00918 { 00919 const System& system = *(pos->second); 00920 const unsigned int nv_sys = system.n_vars(); 00921 00922 system.update_global_solution (sys_soln); 00923 00924 std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element 00925 00926 // Loop over the variable names and load them in order 00927 for (unsigned int var=0; var < nv_sys; ++var) 00928 { 00929 if ( system.variable_type( var ) != type || 00930 ( is_filter_names && std::find(filter_names.begin(), filter_names.end(), system.variable_name( var )) == filter_names.end()) ) 00931 continue; 00932 00933 names.push_back( system.variable_name( var ) ); 00934 00935 // Record the offset for the first element for this variable. 00936 const std::size_t offset = soln.size(); 00937 // Increase size of soln for this variable. 00938 soln.resize( offset + ne ); 00939 00940 const Variable & variable = system.variable(var); 00941 const DofMap & dof_map = system.get_dof_map(); 00942 00943 MeshBase::element_iterator it = _mesh.active_local_elements_begin(); 00944 const MeshBase::element_iterator end_elem = _mesh.active_local_elements_end(); 00945 00946 for ( ; it != end_elem; ++it) 00947 { 00948 if (variable.active_on_subdomain((*it)->subdomain_id())) 00949 { 00950 const Elem* elem = *it; 00951 00952 dof_map.dof_indices (elem, dof_indices, var); 00953 00954 libmesh_assert_equal_to ( 1, dof_indices.size() ); 00955 00956 soln[offset+elem->id()] = sys_soln[dof_indices[0]]; 00957 } 00958 } 00959 00960 } // end loop on variables in this system 00961 00962 } // end loop over systems 00963 00964 // Now each processor has computed contributions to the 00965 // soln vector. Gather them all up. 00966 CommWorld.sum(soln); 00967 }
| System& libMesh::EquationSystems::get_system | ( | const unsigned int | num | ) |
- Returns:
- a writeable referene to the system number
num.
| const System& libMesh::EquationSystems::get_system | ( | const unsigned int | num | ) | const |
- Returns:
- a constant reference to system number
num.
| System& libMesh::EquationSystems::get_system | ( | const std::string & | name | ) |
- Returns:
- a writeable referene to the system named
name.
| const System& libMesh::EquationSystems::get_system | ( | const std::string & | name | ) | const |
- Returns:
- a constant reference to the system named
name.
| System & libMesh::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 587 of file equation_systems.h.
References _systems, end, libMesh::err, and n_systems().
00588 { 00589 libmesh_assert_less (num, this->n_systems()); 00590 00591 const_system_iterator pos = _systems.begin(); 00592 const const_system_iterator end = _systems.end(); 00593 00594 for (; pos != end; ++pos) 00595 if (pos->second->number() == num) 00596 break; 00597 00598 // Check for errors 00599 if (pos == end) 00600 { 00601 libMesh::err << "ERROR: no system number " << num << " found!" 00602 << std::endl; 00603 libmesh_error(); 00604 } 00605 00606 // Attempt dynamic cast 00607 return *libmesh_cast_ptr<T_sys*>(pos->second); 00608 }
| const System & libMesh::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 558 of file equation_systems.h.
References _systems, end, libMesh::err, and n_systems().
00559 { 00560 libmesh_assert_less (num, this->n_systems()); 00561 00562 00563 const_system_iterator pos = _systems.begin(); 00564 const const_system_iterator end = _systems.end(); 00565 00566 for (; pos != end; ++pos) 00567 if (pos->second->number() == num) 00568 break; 00569 00570 // Check for errors 00571 if (pos == end) 00572 { 00573 libMesh::err << "ERROR: no system number " << num << " found!" 00574 << std::endl; 00575 libmesh_error(); 00576 } 00577 00578 // Attempt dynamic cast 00579 return *libmesh_cast_ptr<T_sys*>(pos->second); 00580 }
| System & libMesh::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 640 of file equation_systems.h.
References _systems, and libMesh::err.
00641 { 00642 system_iterator pos = _systems.find(name); 00643 00644 // Check for errors 00645 if (pos == _systems.end()) 00646 { 00647 libMesh::err << "ERROR: no system named " << name << " found!" 00648 << std::endl; 00649 libmesh_error(); 00650 } 00651 00652 // Attempt dynamic cast 00653 return *libmesh_cast_ptr<T_sys*>(pos->second); 00654 }
| const System & libMesh::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 617 of file equation_systems.h.
References _systems, and libMesh::err.
Referenced by libMesh::ExactSolution::_compute_error(), _read_impl(), libMesh::EnsightIO::add_scalar(), add_system(), libMesh::EnsightIO::add_vector(), adjoint_solve(), allgather(), libMesh::ExactSolution::attach_exact_deriv(), libMesh::ExactSolution::attach_exact_hessian(), compare(), libMesh::ExactSolution::compute_error(), libMesh::GMVIO::copy_nodal_solution(), libMesh::ExactSolution::error_norm(), libMesh::ErrorEstimator::estimate_errors(), libMesh::ExactSolution::ExactSolution(), init(), reinit(), sensitivity_solve(), solve(), libMesh::VTKIO::system_vectors_to_vtk(), update(), libMesh::EnsightIO::write_scalar_ascii(), and libMesh::EnsightIO::write_vector_ascii().
00618 { 00619 const_system_iterator pos = _systems.find(name); 00620 00621 // Check for errors 00622 if (pos == _systems.end()) 00623 { 00624 libMesh::err << "ERROR: no system named \"" << name << "\" found!" 00625 << std::endl; 00626 libmesh_error(); 00627 } 00628 00629 // Attempt dynamic cast 00630 return *libmesh_cast_ptr<T_sys*>(pos->second); 00631 }
| bool libMesh::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 495 of file equation_systems.h.
References _mesh_data.
00496 { 00497 return (_mesh_data!=NULL); 00498 }
| bool libMesh::EquationSystems::has_system | ( | const std::string & | name | ) | const [inline] |
- Returns:
- true if the system named
nameexists within this EquationSystems object.
Definition at line 546 of file equation_systems.h.
References _systems.
Referenced by libMesh::EnsightIO::add_scalar(), libMesh::EnsightIO::add_vector(), libMesh::ExactSolution::compute_error(), and libMesh::ExactSolution::error_norm().
| void libMesh::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 163 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
00164 { 00165 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00166 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00167 00168 p.first++; 00169 }
| void libMesh::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 176 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
00177 { 00178 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00179 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00180 00181 p.second++; 00182 }
| void libMesh::EquationSystems::init | ( | ) | [virtual] |
Initialize all the systems
Definition at line 93 of file equation_systems.C.
References _mesh, libMesh::MeshRefinement::clean_refinement_flags(), libMesh::MeshBase::delete_remote_elements(), libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), get_system(), libMesh::n_processors(), n_systems(), libMesh::MeshBase::nodes_begin(), and libMesh::MeshBase::nodes_end().
Referenced by _read_impl(), libMesh::Solver::init(), and libMesh::ErrorVector::plot_error().
00094 { 00095 const unsigned int n_sys = this->n_systems(); 00096 00097 libmesh_assert_not_equal_to (n_sys, 0); 00098 00099 // Distribute the mesh if possible 00100 if (libMesh::n_processors() > 1) 00101 _mesh.delete_remote_elements(); 00102 00103 // Tell all the \p DofObject entities how many systems 00104 // there are. 00105 { 00106 MeshBase::node_iterator node_it = _mesh.nodes_begin(); 00107 const MeshBase::node_iterator node_end = _mesh.nodes_end(); 00108 00109 for ( ; node_it != node_end; ++node_it) 00110 (*node_it)->set_n_systems(n_sys); 00111 00112 MeshBase::element_iterator elem_it = _mesh.elements_begin(); 00113 const MeshBase::element_iterator elem_end = _mesh.elements_end(); 00114 00115 for ( ; elem_it != elem_end; ++elem_it) 00116 (*elem_it)->set_n_systems(n_sys); 00117 } 00118 00119 for (unsigned int i=0; i != this->n_systems(); ++i) 00120 this->get_system(i).init(); 00121 00122 #ifdef LIBMESH_ENABLE_AMR 00123 MeshRefinement mesh_refine(_mesh); 00124 mesh_refine.clean_refinement_flags(); 00125 #endif 00126 }
| std::size_t libMesh::EquationSystems::n_active_dofs | ( | ) | const |
Returns the number of active degrees of freedom for the EquationSystems object.
Definition at line 1218 of file equation_systems.C.
01219 { 01220 std::size_t tot=0; 01221 01222 const_system_iterator pos = _systems.begin(); 01223 const const_system_iterator end = _systems.end(); 01224 01225 for (; pos != end; ++pos) 01226 tot += pos->second->n_active_dofs(); 01227 01228 return tot; 01229 }
| std::size_t libMesh::EquationSystems::n_dofs | ( | ) | const |
- Returns:
- the total number of degrees of freedom in all systems.
Definition at line 1202 of file equation_systems.C.
01203 { 01204 std::size_t tot=0; 01205 01206 const_system_iterator pos = _systems.begin(); 01207 const const_system_iterator end = _systems.end(); 01208 01209 for (; pos != end; ++pos) 01210 tot += pos->second->n_dofs(); 01211 01212 return tot; 01213 }
| static unsigned int libMesh::ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 79 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
00080 { return _n_objects; }
| unsigned int libMesh::EquationSystems::n_systems | ( | ) | const [inline] |
- Returns:
- the number of equation systems.
Definition at line 502 of file equation_systems.h.
References _systems.
Referenced by add_system(), adjoint_solve(), allgather(), libMesh::ExactSolution::attach_exact_deriv(), libMesh::ExactSolution::attach_exact_hessian(), build_discontinuous_solution_vector(), build_solution_vector(), build_variable_names(), compare(), libMesh::GMVIO::copy_nodal_solution(), libMesh::ErrorEstimator::estimate_errors(), libMesh::ExactSolution::ExactSolution(), get_info(), get_solution(), get_system(), init(), reinit(), sensitivity_solve(), solve(), libMesh::VTKIO::system_vectors_to_vtk(), update(), and write().
00503 { 00504 return libmesh_cast_int<unsigned int>(_systems.size()); 00505 }
| unsigned int libMesh::EquationSystems::n_vars | ( | ) | const |
- Returns:
- the total number of variables in all systems.
Definition at line 1187 of file equation_systems.C.
Referenced by build_discontinuous_solution_vector(), and build_variable_names().
01188 { 01189 unsigned int tot=0; 01190 01191 const_system_iterator pos = _systems.begin(); 01192 const const_system_iterator end = _systems.end(); 01193 01194 for (; pos != end; ++pos) 01195 tot += pos->second->n_vars(); 01196 01197 return tot; 01198 }
| void libMesh::ReferenceCounter::print_info | ( | std::ostream & | out = libMesh::out |
) | [static, inherited] |
Prints the reference information, by default to libMesh::out.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
00089 { 00090 if( _enable_print_counter ) out_stream << ReferenceCounter::get_info(); 00091 }
| void libMesh::EquationSystems::print_info | ( | std::ostream & | os = libMesh::out |
) | const |
Prints information about the equation systems, by default to libMesh::out.
Definition at line 1171 of file equation_systems.C.
References get_info().
Referenced by libMesh::operator<<().
01172 { 01173 os << this->get_info() 01174 << std::endl; 01175 }
| void libMesh::EquationSystems::read | ( | const std::string & | name, | |
| const unsigned int | read_flags = (READ_HEADER | READ_DATA) | |||
| ) |
Definition at line 72 of file equation_systems_io.C.
References _mesh, libMesh::MeshRefinement::clean_refinement_flags(), libMeshEnums::DECODE, read(), and libMeshEnums::READ.
00074 { 00075 libMeshEnums::XdrMODE mode = READ; 00076 if (name.find(".xdr") != std::string::npos) 00077 mode = DECODE; 00078 this->read(name, mode, read_flags); 00079 00080 #ifdef LIBMESH_ENABLE_AMR 00081 MeshRefinement mesh_refine(_mesh); 00082 mesh_refine.clean_refinement_flags(); 00083 #endif 00084 }
| void libMesh::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 88 of file equation_systems_io.C.
References _mesh, _read_impl(), libMesh::MeshRefinement::clean_refinement_flags(), libMesh::out, and TRY_READ_IFEMS.
Referenced by _read_impl(), and read().
00091 { 00092 #ifdef LIBMESH_ENABLE_EXCEPTIONS 00093 00094 // If we have exceptions enabled we can be considerate and try 00095 // to read old restart files which contain infinite element 00096 // information but do not have the " with infinite elements" 00097 // string in the version information. 00098 00099 // First try the read the user requested 00100 try 00101 { 00102 this->_read_impl (name, mode, read_flags); 00103 } 00104 00105 // If that fails, try it again but explicitly request we look for infinite element info 00106 catch (...) 00107 { 00108 libMesh::out << "\n*********************************************************************\n" 00109 << "READING THE FILE \"" << name << "\" FAILED.\n" 00110 << "It is possible this file contains infinite element information,\n" 00111 << "but the version string does not contain \" with infinite elements\"\n" 00112 << "Let's try this again, but looking for infinite element information...\n" 00113 << "*********************************************************************\n" 00114 << std::endl; 00115 00116 try 00117 { 00118 this->_read_impl (name, mode, read_flags | EquationSystems::TRY_READ_IFEMS); 00119 } 00120 00121 // If all that failed, we are out of ideas here... 00122 catch (...) 00123 { 00124 libMesh::out << "\n*********************************************************************\n" 00125 << "Well, at least we tried!\n" 00126 << "Good Luck!!\n" 00127 << "*********************************************************************\n" 00128 << std::endl; 00129 throw; 00130 } 00131 } 00132 00133 #else 00134 00135 // no exceptions - cross your fingers... 00136 this->_read_impl (name, mode, read_flags); 00137 00138 #endif // #ifdef LIBMESH_ENABLE_EXCEPTIONS 00139 00140 #ifdef LIBMESH_ENABLE_AMR 00141 MeshRefinement mesh_refine(_mesh); 00142 mesh_refine.clean_refinement_flags(); 00143 #endif 00144 }
| void libMesh::EquationSystems::reinit | ( | ) | [virtual] |
Reinitialize all the systems
Definition at line 130 of file equation_systems.C.
References _mesh, libMesh::MeshRefinement::coarsen_elements(), libMesh::MeshBase::contract(), libMesh::DofMap::create_dof_constraints(), libMesh::DofMap::distribute_dofs(), libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), libMesh::MeshRefinement::face_level_mismatch_limit(), libMesh::System::get_dof_map(), get_mesh(), get_system(), libMesh::DofObject::n_systems(), n_systems(), libMesh::System::n_vars(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), libMesh::DofMap::prepare_send_list(), libMesh::DofMap::process_constraints(), libMesh::System::prolong_vectors(), libMesh::MeshRefinement::refine_elements(), libMesh::System::restrict_vectors(), libMesh::System::time, and libMesh::System::user_constrain().
Referenced by libMesh::AdjointRefinementEstimator::estimate_error().
00131 { 00132 parallel_only(); 00133 00134 libmesh_assert_not_equal_to (this->n_systems(), 0); 00135 00136 #ifdef DEBUG 00137 // Make sure all the \p DofObject entities know how many systems 00138 // there are. 00139 { 00140 // All the nodes 00141 MeshBase::node_iterator node_it = _mesh.nodes_begin(); 00142 const MeshBase::node_iterator node_end = _mesh.nodes_end(); 00143 00144 for ( ; node_it != node_end; ++node_it) 00145 { 00146 Node *node = *node_it; 00147 libmesh_assert_equal_to (node->n_systems(), this->n_systems()); 00148 } 00149 00150 // All the elements 00151 MeshBase::element_iterator elem_it = _mesh.elements_begin(); 00152 const MeshBase::element_iterator elem_end = _mesh.elements_end(); 00153 00154 for ( ; elem_it != elem_end; ++elem_it) 00155 { 00156 Elem *elem = *elem_it; 00157 libmesh_assert_equal_to (elem->n_systems(), this->n_systems()); 00158 } 00159 } 00160 #endif 00161 00162 // Localize each system's vectors 00163 for (unsigned int i=0; i != this->n_systems(); ++i) 00164 this->get_system(i).re_update(); 00165 00166 #ifdef LIBMESH_ENABLE_AMR 00167 00168 bool dof_constraints_created = false; 00169 bool mesh_changed = false; 00170 00171 // FIXME: For backwards compatibility, assume 00172 // refine_and_coarsen_elements or refine_uniformly have already 00173 // been called 00174 { 00175 for (unsigned int i=0; i != this->n_systems(); ++i) 00176 { 00177 System &sys = this->get_system(i); 00178 00179 // Don't do anything if the system doesn't have any variables in it 00180 if(!sys.n_vars()) 00181 continue; 00182 00183 sys.get_dof_map().distribute_dofs(_mesh); 00184 00185 // Recreate any hanging node constraints 00186 sys.get_dof_map().create_dof_constraints(_mesh, sys.time); 00187 00188 // Apply any user-defined constraints 00189 sys.user_constrain(); 00190 00191 // Expand any recursive constraints 00192 sys.get_dof_map().process_constraints(_mesh); 00193 00194 // And clean up the send_list before we use it again 00195 sys.get_dof_map().prepare_send_list(); 00196 00197 sys.prolong_vectors(); 00198 } 00199 mesh_changed = true; 00200 dof_constraints_created = true; 00201 } 00202 00203 // FIXME: Where should the user set maintain_level_one now?? 00204 // Don't override previous settings, for now 00205 00206 MeshRefinement mesh_refine(_mesh); 00207 00208 mesh_refine.face_level_mismatch_limit() = false; 00209 00210 // Try to coarsen the mesh, then restrict each system's vectors 00211 // if necessary 00212 if (mesh_refine.coarsen_elements()) 00213 { 00214 for (unsigned int i=0; i != this->n_systems(); ++i) 00215 { 00216 System &sys = this->get_system(i); 00217 if (!dof_constraints_created) 00218 { 00219 sys.get_dof_map().distribute_dofs(_mesh); 00220 sys.get_dof_map().create_dof_constraints(_mesh, sys.time); 00221 sys.user_constrain(); 00222 sys.get_dof_map().process_constraints(_mesh); 00223 sys.get_dof_map().prepare_send_list(); 00224 00225 } 00226 sys.restrict_vectors(); 00227 } 00228 mesh_changed = true; 00229 dof_constraints_created = true; 00230 } 00231 00232 // Once vectors are all restricted, we can delete 00233 // children of coarsened elements 00234 if (mesh_changed) 00235 this->get_mesh().contract(); 00236 00237 // Try to refine the mesh, then prolong each system's vectors 00238 // if necessary 00239 if (mesh_refine.refine_elements()) 00240 { 00241 for (unsigned int i=0; i != this->n_systems(); ++i) 00242 { 00243 System &sys = this->get_system(i); 00244 if (!dof_constraints_created) 00245 { 00246 sys.get_dof_map().distribute_dofs(_mesh); 00247 sys.get_dof_map().create_dof_constraints(_mesh, sys.time); 00248 sys.user_constrain(); 00249 sys.get_dof_map().process_constraints(_mesh); 00250 sys.get_dof_map().prepare_send_list(); 00251 00252 } 00253 sys.prolong_vectors(); 00254 } 00255 mesh_changed = true; 00256 dof_constraints_created = true; 00257 } 00258 00259 // If the mesh has changed, systems will need to create new dof 00260 // constraints and update their global solution vectors 00261 if (mesh_changed) 00262 { 00263 for (unsigned int i=0; i != this->n_systems(); ++i) 00264 this->get_system(i).reinit(); 00265 } 00266 #endif // #ifdef LIBMESH_ENABLE_AMR 00267 }
| void libMesh::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 452 of file equation_systems.C.
References get_system(), and n_systems().
00453 { 00454 libmesh_assert (this->n_systems()); 00455 00456 for (unsigned int i=0; i != this->n_systems(); ++i) 00457 this->get_system(i).sensitivity_solve(parameters_in); 00458 }
| void libMesh::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 442 of file equation_systems.C.
References get_system(), and n_systems().
Referenced by libMesh::Solver::solve().
00443 { 00444 libmesh_assert (this->n_systems()); 00445 00446 for (unsigned int i=0; i != this->n_systems(); ++i) 00447 this->get_system(i).solve(); 00448 }
| void libMesh::EquationSystems::update | ( | ) |
Updates local values for all the systems
Definition at line 322 of file equation_systems.C.
References get_system(), and n_systems().
Referenced by _read_impl().
00323 { 00324 START_LOG("update()","EquationSystems"); 00325 00326 // Localize each system's vectors 00327 for (unsigned int i=0; i != this->n_systems(); ++i) 00328 this->get_system(i).update(); 00329 00330 STOP_LOG("update()","EquationSystems"); 00331 }
| void libMesh::EquationSystems::write | ( | const std::string & | name, | |
| const unsigned int | write_flags = (WRITE_DATA) | |||
| ) | const |
Definition at line 370 of file equation_systems_io.C.
References libMeshEnums::ENCODE, write(), and libMeshEnums::WRITE.
| void libMesh::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 381 of file equation_systems_io.C.
References _mesh, _systems, libMesh::Xdr::data(), libMesh::get_io_compatibility_version(), get_mesh(), libMesh::MeshTools::Private::globally_renumber_nodes_and_elements(), libMesh::MeshBase::is_serial(), mesh, n_systems(), libMesh::processor_id(), libMesh::Xdr::set_version(), WRITE_ADDITIONAL_DATA, WRITE_DATA, WRITE_PARALLEL_FILES, WRITE_SERIAL_FILES, and libMesh::Xdr::writing().
Referenced by write().
00384 { 00448 // the EquationSystems::write() method should look constant, 00449 // but we need to assign a temporary numbering to the nodes 00450 // and elements in the mesh, which requires that we abuse const_cast 00451 { 00452 MeshBase &mesh = const_cast<MeshBase&>(this->get_mesh()); 00453 MeshTools::Private::globally_renumber_nodes_and_elements(mesh); 00454 } 00455 00456 // set booleans from write_flags argument 00457 const bool write_data = write_flags & EquationSystems::WRITE_DATA; 00458 const bool write_additional_data = write_flags & EquationSystems::WRITE_ADDITIONAL_DATA; 00459 00460 // always write parallel files if we're instructed to write in 00461 // parallel 00462 const bool write_parallel_files = 00463 (write_flags & EquationSystems::WRITE_PARALLEL_FILES) || 00464 // but also write parallel files if we haven't been instructed to 00465 // write in serial and we're on a distributed mesh 00466 (!(write_flags & EquationSystems::WRITE_SERIAL_FILES) && 00467 !this->get_mesh().is_serial()); 00468 00469 // New scope so that io will close before we try to zip the file 00470 { 00471 Xdr io((libMesh::processor_id()==0) ? name : "", mode); 00472 libmesh_assert (io.writing()); 00473 00474 START_LOG("write()","EquationSystems"); 00475 00476 const unsigned int proc_id = libMesh::processor_id(); 00477 unsigned int n_sys = this->n_systems(); 00478 00479 std::map<std::string, System*>::const_iterator 00480 pos = _systems.begin(); 00481 00482 std::string comment; 00483 char buf[256]; 00484 00485 // set the version number in the Xdr object 00486 io.set_version(LIBMESH_VERSION_ID(LIBMESH_MAJOR_VERSION, 00487 LIBMESH_MINOR_VERSION, 00488 LIBMESH_MICRO_VERSION)); 00489 00490 // Only write the header information 00491 // if we are processor 0. 00492 if (proc_id == 0) 00493 { 00494 // 1.) 00495 // Write the version header 00496 std::string version("libMesh-" + libMesh::get_io_compatibility_version()); 00497 if (write_parallel_files) version += " parallel"; 00498 00499 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00500 version += " with infinite elements"; 00501 #endif 00502 io.data (version, "# File Format Identifier"); 00503 00504 // 2.) 00505 // Write the number of equation systems 00506 io.data (n_sys, "# No. of Equation Systems"); 00507 00508 while (pos != _systems.end()) 00509 { 00510 // 3.) 00511 // Write the name of the sys_num-th system 00512 { 00513 const unsigned int sys_num = pos->second->number(); 00514 std::string sys_name = pos->first; 00515 00516 comment = "# Name, System No. "; 00517 std::sprintf(buf, "%d", sys_num); 00518 comment += buf; 00519 00520 io.data (sys_name, comment.c_str()); 00521 } 00522 00523 // 4.) 00524 // Write the type of system handled 00525 { 00526 const unsigned int sys_num = pos->second->number(); 00527 std::string sys_type = pos->second->system_type(); 00528 00529 comment = "# Type, System No. "; 00530 std::sprintf(buf, "%d", sys_num); 00531 comment += buf; 00532 00533 io.data (sys_type, comment.c_str()); 00534 } 00535 00536 // 5.) - 9.) 00537 // Let System::write_header() do the job 00538 pos->second->write_header (io, version, write_additional_data); 00539 00540 ++pos; 00541 } 00542 } 00543 00544 // Start from the first system, again, 00545 // to write vectors to disk, if wanted 00546 if (write_data) 00547 { 00548 // open a parallel buffer if warranted. 00549 Xdr local_io (write_parallel_files ? local_file_name(name) : "", mode); 00550 00551 for (pos = _systems.begin(); pos != _systems.end(); ++pos) 00552 { 00553 // 10.) + 11.) 00554 if (write_parallel_files) 00555 pos->second->write_parallel_data (local_io,write_additional_data); 00556 else 00557 pos->second->write_serialized_data (io,write_additional_data); 00558 } 00559 } 00560 00561 STOP_LOG("write()","EquationSystems"); 00562 } 00563 00564 // the EquationSystems::write() method should look constant, 00565 // but we need to undo the temporary numbering of the nodes 00566 // and elements in the mesh, which requires that we abuse const_cast 00567 const_cast<MeshBase&>(_mesh).fix_broken_node_and_element_numbering(); 00568 }
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.
Member Data Documentation
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited] |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 137 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
MeshBase& libMesh::EquationSystems::_mesh [protected] |
The mesh data structure
Definition at line 418 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(), get_solution(), init(), read(), reinit(), and write().
MeshData* libMesh::EquationSystems::_mesh_data [protected] |
A pointer to the MeshData object you would like to use. Can be NULL.
Definition at line 424 of file equation_systems.h.
Referenced by get_mesh_data(), and has_mesh_data().
Threads::spin_mutex libMesh::ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 131 of file reference_counter.h.
Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 126 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
std::map<std::string, System*> libMesh::EquationSystems::_systems [protected] |
Data structure holding the systems.
Definition at line 429 of file equation_systems.h.
Referenced by add_system(), build_discontinuous_solution_vector(), build_solution_vector(), build_variable_names(), clear(), compare(), delete_system(), get_info(), get_solution(), get_system(), has_system(), n_active_dofs(), n_dofs(), n_systems(), n_vars(), and write().
Data structure holding arbitrary parameters.
Definition at line 409 of file equation_systems.h.
Referenced by libMesh::ExactSolution::attach_exact_deriv(), libMesh::ExactSolution::attach_exact_hessian(), libMesh::NewmarkSystem::clear(), clear(), libMesh::FrequencySystem::clear_all(), EquationSystems(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::FEComputeData::init(), libMesh::FrequencySystem::init_data(), libMesh::FrequencySystem::n_frequencies(), libMesh::NewmarkSystem::NewmarkSystem(), libMesh::NonlinearImplicitSystem::NonlinearImplicitSystem(), libMesh::FrequencySystem::set_current_frequency(), libMesh::FrequencySystem::set_frequencies(), libMesh::FrequencySystem::set_frequencies_by_range(), libMesh::FrequencySystem::set_frequencies_by_steps(), libMesh::NewmarkSystem::set_newmark_parameters(), libMesh::NonlinearImplicitSystem::set_solver_parameters(), libMesh::LinearImplicitSystem::solve(), libMesh::FrequencySystem::solve(), libMesh::EigenSystem::solve(), libMesh::CondensedEigenSystem::solve(), and libMesh::WrappedFunction< Output >::WrappedFunction().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:17 UTC
Hosted By: