$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 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 00020 // C++ includes 00021 00022 00023 // Local includes 00024 #include "libmesh/system_subset_by_subdomain.h" 00025 #include "libmesh/system.h" 00026 #include "libmesh/dof_map.h" 00027 #include "libmesh/parallel.h" 00028 00029 namespace libMesh 00030 { 00031 // ------------------------------------------------------------ 00032 // SubdomainSelection implementation 00033 SystemSubsetBySubdomain::SubdomainSelection:: 00034 SubdomainSelection (void) 00035 { 00036 } 00037 00038 SystemSubsetBySubdomain::SubdomainSelection:: 00039 ~SubdomainSelection (void) 00040 { 00041 } 00042 00043 SystemSubsetBySubdomain::SubdomainSelectionByList:: 00044 SubdomainSelectionByList (const std::set<subdomain_id_type>& list): 00045 _list(list) 00046 { 00047 } 00048 00049 bool 00050 SystemSubsetBySubdomain::SubdomainSelectionByList:: 00051 operator()(const subdomain_id_type& subdomain_id)const 00052 { 00053 return _list.find(subdomain_id)!=_list.end(); 00054 } 00055 00056 // ------------------------------------------------------------ 00057 // SystemSubsetBySubdomain implementation 00058 00059 SystemSubsetBySubdomain:: 00060 SystemSubsetBySubdomain (const System& system, 00061 const SubdomainSelection& subdomain_selection, 00062 const std::set<unsigned int>* const var_nums): 00063 SystemSubset(system), 00064 ParallelObject(system), 00065 _var_nums(), 00066 _dof_ids() 00067 { 00068 this->set_var_nums(var_nums); 00069 this->init(subdomain_selection); 00070 } 00071 00072 SystemSubsetBySubdomain:: 00073 SystemSubsetBySubdomain (const System& system, 00074 const std::set<subdomain_id_type>& subdomain_ids, 00075 const std::set<unsigned int>* const var_nums): 00076 SystemSubset(system), 00077 ParallelObject(system), 00078 _var_nums(), 00079 _dof_ids() 00080 { 00081 this->set_var_nums(var_nums); 00082 this->init(subdomain_ids); 00083 } 00084 00085 SystemSubsetBySubdomain:: 00086 ~SystemSubsetBySubdomain (void) 00087 { 00088 } 00089 00090 const std::vector<unsigned int>& 00091 SystemSubsetBySubdomain:: 00092 dof_ids(void)const 00093 { 00094 return _dof_ids; 00095 } 00096 00097 void 00098 SystemSubsetBySubdomain:: 00099 set_var_nums (const std::set<unsigned int>* const var_nums) 00100 { 00101 _var_nums.clear(); 00102 if(var_nums!=NULL) 00103 { 00104 _var_nums = *var_nums; 00105 } 00106 else 00107 { 00108 for(unsigned int i=0; i<_system.n_vars(); i++) 00109 { 00110 _var_nums.insert(i); 00111 } 00112 } 00113 } 00114 00115 void 00116 SystemSubsetBySubdomain:: 00117 init (const SubdomainSelection& subdomain_selection) 00118 { 00119 _dof_ids.clear(); 00120 00121 std::vector<std::vector<dof_id_type> > dof_ids_per_processor(this->n_processors()); 00122 00123 const DofMap & dof_map = _system.get_dof_map(); 00124 std::vector<dof_id_type> dof_indices; 00125 00126 const MeshBase& mesh = _system.get_mesh(); 00127 MeshBase::const_element_iterator el = mesh.active_local_elements_begin(); 00128 const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end(); 00129 for ( ; el != end_el; ++el) 00130 { 00131 const Elem* elem = *el; 00132 if(subdomain_selection(elem->subdomain_id())) 00133 { 00134 std::set<unsigned int>::const_iterator it = _var_nums.begin(); 00135 const std::set<unsigned int>::const_iterator itEnd = _var_nums.end(); 00136 for (; it!=itEnd; ++it) 00137 { 00138 dof_map.dof_indices (elem, dof_indices, *it); 00139 for(size_t i=0; i<dof_indices.size(); i++) 00140 { 00141 const dof_id_type dof = dof_indices[i]; 00142 for(processor_id_type proc=0; proc<this->n_processors(); proc++) 00143 { 00144 if((dof>=dof_map.first_dof(proc)) && (dof<dof_map.end_dof(proc))) 00145 { 00146 dof_ids_per_processor[proc].push_back(dof); 00147 } 00148 } 00149 } 00150 } 00151 } 00152 } 00153 00154 /* Distribute information among processors. */ 00155 std::vector<Parallel::Request> request_per_processor(this->n_processors()); 00156 for(unsigned int proc=0; proc<this->n_processors(); proc++) 00157 { 00158 if(proc!=this->processor_id()) 00159 { 00160 this->comm().send(proc,dof_ids_per_processor[proc],request_per_processor[proc]); 00161 } 00162 } 00163 for(unsigned int proc=0; proc<this->n_processors(); proc++) 00164 { 00165 std::vector<dof_id_type> received_dofs; 00166 if(proc==this->processor_id()) 00167 { 00168 received_dofs = dof_ids_per_processor[proc]; 00169 } 00170 else 00171 { 00172 this->comm().receive(proc,received_dofs); 00173 } 00174 for(unsigned int i=0; i<received_dofs.size(); i++) 00175 { 00176 _dof_ids.push_back(received_dofs[i]); 00177 } 00178 } 00179 00180 /* Sort and unique the vector (using the same mechanism as in \p 00181 DofMap::prepare_send_list()). */ 00182 std::sort(_dof_ids.begin(), _dof_ids.end()); 00183 std::vector<unsigned int>::iterator new_end = std::unique (_dof_ids.begin(), _dof_ids.end()); 00184 std::vector<unsigned int> (_dof_ids.begin(), new_end).swap (_dof_ids); 00185 00186 /* Wait for sends to be complete. */ 00187 for(unsigned int proc=0; proc<this->n_processors(); proc++) 00188 { 00189 if(proc!=this->processor_id()) 00190 { 00191 request_per_processor[proc].wait(); 00192 } 00193 } 00194 } 00195 00196 void 00197 SystemSubsetBySubdomain:: 00198 init (const std::set<subdomain_id_type>& subdomain_ids) 00199 { 00200 SubdomainSelectionByList selection(subdomain_ids); 00201 this->init(selection); 00202 } 00203 00204 } // namespace libMesh