equation_systems.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 // System includes 00020 #include <sstream> 00021 00022 // Local Includes 00023 #include "libmesh/explicit_system.h" 00024 #include "libmesh/fe_interface.h" 00025 #include "libmesh/frequency_system.h" 00026 #include "libmesh/linear_implicit_system.h" 00027 #include "libmesh/mesh_refinement.h" 00028 #include "libmesh/newmark_system.h" 00029 #include "libmesh/nonlinear_implicit_system.h" 00030 #include "libmesh/rb_construction.h" 00031 #include "libmesh/transient_rb_construction.h" 00032 #include "libmesh/eigen_system.h" 00033 #include "libmesh/parallel.h" 00034 #include "libmesh/transient_system.h" 00035 #include "libmesh/dof_map.h" 00036 #include "libmesh/mesh_base.h" 00037 #include "libmesh/elem.h" 00038 #include "libmesh/libmesh_logging.h" 00039 00040 // Include the systems before this one to avoid 00041 // overlapping forward declarations. 00042 #include "libmesh/equation_systems.h" 00043 00044 namespace libMesh 00045 { 00046 00047 // Forward Declarations 00048 00049 00050 00051 00052 // ------------------------------------------------------------ 00053 // EquationSystems class implementation 00054 EquationSystems::EquationSystems (MeshBase& m, MeshData* mesh_data) : 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 } 00062 00063 00064 00065 EquationSystems::~EquationSystems () 00066 { 00067 this->clear (); 00068 } 00069 00070 00071 00072 void EquationSystems::clear () 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 } 00090 00091 00092 00093 void EquationSystems::init () 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 } 00127 00128 00129 00130 void EquationSystems::reinit () 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 } 00268 00269 00270 00271 void EquationSystems::allgather () 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 } 00318 00319 00320 00321 00322 void EquationSystems::update () 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 } 00332 00333 00334 00335 System & EquationSystems::add_system (const std::string& sys_type, 00336 const std::string& name) 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 } 00417 00418 00419 00420 00421 00422 00423 void EquationSystems::delete_system (const std::string& name) 00424 { 00425 libmesh_deprecated(); 00426 00427 if (!_systems.count(name)) 00428 { 00429 libMesh::err << "ERROR: no system named " 00430 << name << std::endl; 00431 00432 libmesh_error(); 00433 } 00434 00435 delete _systems[name]; 00436 00437 _systems.erase (name); 00438 } 00439 00440 00441 00442 void EquationSystems::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 } 00449 00450 00451 00452 void EquationSystems::sensitivity_solve (const ParameterVector& parameters_in) 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 } 00459 00460 00461 00462 void EquationSystems::adjoint_solve (const QoISet& qoi_indices) 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 } 00469 00470 00471 00472 void EquationSystems::build_variable_names (std::vector<std::string>& var_names, 00473 const FEType *type, 00474 const std::set<std::string>* system_names) const 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 } 00583 00584 00585 00586 void EquationSystems::build_solution_vector (std::vector<Number>&, 00587 const std::string&, 00588 const std::string&) const 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 } 00669 00670 00671 00672 00673 void EquationSystems::build_solution_vector (std::vector<Number>& soln, 00674 const std::set<std::string>* system_names) const 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 } 00883 00884 00885 void EquationSystems::get_solution (std::vector<Number>& soln, 00886 std::vector<std::string> & names ) const 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 } 00968 00969 00970 00971 void EquationSystems::build_discontinuous_solution_vector (std::vector<Number>& soln) const 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 } 01069 01070 01071 01072 bool EquationSystems::compare (const EquationSystems& other_es, 01073 const Real threshold, 01074 const bool verbose) const 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 } 01130 01131 01132 01133 std::string EquationSystems::get_info () const 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 } 01168 01169 01170 01171 void EquationSystems::print_info (std::ostream& os) const 01172 { 01173 os << this->get_info() 01174 << std::endl; 01175 } 01176 01177 01178 01179 std::ostream& operator << (std::ostream& os, const EquationSystems& es) 01180 { 01181 es.print_info(os); 01182 return os; 01183 } 01184 01185 01186 01187 unsigned int EquationSystems::n_vars () const 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 } 01199 01200 01201 01202 std::size_t EquationSystems::n_dofs () const 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 } 01214 01215 01216 01217 01218 std::size_t EquationSystems::n_active_dofs () const 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 } 01230 01231 01232 void EquationSystems::_add_system_to_nodes_and_elems() 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 } 01248 01249 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: