$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 #include <set> 00022 #include <algorithm> // for std::fill, std::equal_range, std::max, std::lower_bound, etc. 00023 #include <sstream> 00024 00025 // Local Includes ----------------------------------- 00026 #include "libmesh/coupling_matrix.h" 00027 #include "libmesh/dense_matrix.h" 00028 #include "libmesh/dense_vector_base.h" 00029 #include "libmesh/dirichlet_boundaries.h" 00030 #include "libmesh/dof_map.h" 00031 #include "libmesh/elem.h" 00032 #include "libmesh/fe_interface.h" 00033 #include "libmesh/fe_type.h" 00034 #include "libmesh/fe_base.h" // FEBase::build() for continuity test 00035 #include "libmesh/libmesh_logging.h" 00036 #include "libmesh/mesh_base.h" 00037 #include "libmesh/mesh_tools.h" 00038 #include "libmesh/numeric_vector.h" 00039 #include "libmesh/parallel.h" 00040 #include "libmesh/periodic_boundaries.h" 00041 #include "libmesh/sparse_matrix.h" 00042 #include "libmesh/sparsity_pattern.h" 00043 #include "libmesh/string_to_enum.h" 00044 #include "libmesh/threads.h" 00045 #include "libmesh/mesh_subdivision_support.h" 00046 00047 00048 00049 namespace libMesh 00050 { 00051 00052 // ------------------------------------------------------------ 00053 // DofMap member functions 00054 UniquePtr<SparsityPattern::Build> DofMap::build_sparsity 00055 (const MeshBase& mesh) const 00056 { 00057 libmesh_assert (mesh.is_prepared()); 00058 00059 START_LOG("build_sparsity()", "DofMap"); 00060 00061 // Compute the sparsity structure of the global matrix. This can be 00062 // fed into a PetscMatrix to allocate exacly the number of nonzeros 00063 // necessary to store the matrix. This algorithm should be linear 00064 // in the (# of elements)*(# nodes per element) 00065 00066 // We can be more efficient in the threaded sparsity pattern assembly 00067 // if we don't need the exact pattern. For some sparse matrix formats 00068 // a good upper bound will suffice. 00069 00070 // See if we need to include sparsity pattern entries for coupling 00071 // between neighbor dofs 00072 bool implicit_neighbor_dofs = this->use_coupled_neighbor_dofs(mesh); 00073 00074 // We can compute the sparsity pattern in parallel on multiple 00075 // threads. The goal is for each thread to compute the full sparsity 00076 // pattern for a subset of elements. These sparsity patterns can 00077 // be efficiently merged in the SparsityPattern::Build::join() 00078 // method, especially if there is not too much overlap between them. 00079 // Even better, if the full sparsity pattern is not needed then 00080 // the number of nonzeros per row can be estimated from the 00081 // sparsity patterns created on each thread. 00082 UniquePtr<SparsityPattern::Build> sp 00083 (new SparsityPattern::Build (mesh, 00084 *this, 00085 this->_dof_coupling, 00086 implicit_neighbor_dofs, 00087 need_full_sparsity_pattern)); 00088 00089 Threads::parallel_reduce (ConstElemRange (mesh.active_local_elements_begin(), 00090 mesh.active_local_elements_end()), *sp); 00091 00092 sp->parallel_sync(); 00093 00094 #ifndef NDEBUG 00095 // Avoid declaring these variables unless asserts are enabled. 00096 const processor_id_type proc_id = mesh.processor_id(); 00097 const dof_id_type n_dofs_on_proc = this->n_dofs_on_processor(proc_id); 00098 #endif 00099 libmesh_assert_equal_to (sp->sparsity_pattern.size(), n_dofs_on_proc); 00100 00101 STOP_LOG("build_sparsity()", "DofMap"); 00102 00103 // Check to see if we have any extra stuff to add to the sparsity_pattern 00104 if (_extra_sparsity_function) 00105 { 00106 if (_augment_sparsity_pattern) 00107 { 00108 libmesh_here(); 00109 libMesh::out << "WARNING: You have specified both an extra sparsity function and object.\n" 00110 << " Are you sure this is what you meant to do??" 00111 << std::endl; 00112 } 00113 00114 _extra_sparsity_function 00115 (sp->sparsity_pattern, sp->n_nz, 00116 sp->n_oz, _extra_sparsity_context); 00117 } 00118 00119 if (_augment_sparsity_pattern) 00120 _augment_sparsity_pattern->augment_sparsity_pattern 00121 (sp->sparsity_pattern, sp->n_nz, sp->n_oz); 00122 00123 return UniquePtr<SparsityPattern::Build>(sp.release()); 00124 } 00125 00126 00127 00128 DofMap::DofMap(const unsigned int number, 00129 const ParallelObject &parent_decomp) : 00130 ParallelObject (parent_decomp), 00131 _dof_coupling(NULL), 00132 _variables(), 00133 _variable_groups(), 00134 _sys_number(number), 00135 _matrices(), 00136 _first_df(), 00137 _end_df(), 00138 _first_scalar_df(), 00139 _send_list(), 00140 _augment_sparsity_pattern(NULL), 00141 _extra_sparsity_function(NULL), 00142 _extra_sparsity_context(NULL), 00143 _augment_send_list(NULL), 00144 _extra_send_list_function(NULL), 00145 _extra_send_list_context(NULL), 00146 need_full_sparsity_pattern(false), 00147 _n_nz(NULL), 00148 _n_oz(NULL), 00149 _n_dfs(0), 00150 _n_SCALAR_dofs(0) 00151 #ifdef LIBMESH_ENABLE_AMR 00152 , _n_old_dfs(0), 00153 _first_old_df(), 00154 _end_old_df(), 00155 _first_old_scalar_df() 00156 #endif 00157 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00158 , _dof_constraints() 00159 , _stashed_dof_constraints() 00160 , _primal_constraint_values() 00161 , _adjoint_constraint_values() 00162 #endif 00163 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 00164 , _node_constraints() 00165 #endif 00166 #ifdef LIBMESH_ENABLE_PERIODIC 00167 , _periodic_boundaries(new PeriodicBoundaries) 00168 #endif 00169 #ifdef LIBMESH_ENABLE_DIRICHLET 00170 , _dirichlet_boundaries(new DirichletBoundaries) 00171 , _adjoint_dirichlet_boundaries() 00172 #endif 00173 , _implicit_neighbor_dofs_initialized(false), 00174 _implicit_neighbor_dofs(false) 00175 { 00176 _matrices.clear(); 00177 } 00178 00179 00180 00181 // Destructor 00182 DofMap::~DofMap() 00183 { 00184 this->clear(); 00185 #ifdef LIBMESH_ENABLE_PERIODIC 00186 delete _periodic_boundaries; 00187 #endif 00188 #ifdef LIBMESH_ENABLE_DIRICHLET 00189 delete _dirichlet_boundaries; 00190 for (unsigned int q = 0; q != _adjoint_dirichlet_boundaries.size(); ++q) 00191 delete _adjoint_dirichlet_boundaries[q]; 00192 #endif 00193 } 00194 00195 00196 #ifdef LIBMESH_ENABLE_PERIODIC 00197 00198 bool DofMap::is_periodic_boundary (const boundary_id_type boundaryid) const 00199 { 00200 if (_periodic_boundaries->count(boundaryid) != 0) 00201 return true; 00202 00203 return false; 00204 } 00205 00206 #endif 00207 00208 00209 00210 // void DofMap::add_variable (const Variable &var) 00211 // { 00212 // libmesh_not_implemented(); 00213 // _variables.push_back (var); 00214 // } 00215 00216 00217 00218 void DofMap::add_variable_group (const VariableGroup &var_group) 00219 { 00220 _variable_groups.push_back(var_group); 00221 00222 VariableGroup &new_var_group = _variable_groups.back(); 00223 00224 for (unsigned int var=0; var<new_var_group.n_variables(); var++) 00225 _variables.push_back (new_var_group(var)); 00226 } 00227 00228 00229 00230 void DofMap::attach_matrix (SparseMatrix<Number>& matrix) 00231 { 00232 parallel_object_only(); 00233 00234 // We shouldn't be trying to re-attach the same matrices repeatedly 00235 libmesh_assert (std::find(_matrices.begin(), _matrices.end(), 00236 &matrix) == _matrices.end()); 00237 00238 _matrices.push_back(&matrix); 00239 00240 matrix.attach_dof_map (*this); 00241 00242 // If we've already computed sparsity, then it's too late 00243 // to wait for "compute_sparsity" to help with sparse matrix 00244 // initialization, and we need to handle this matrix individually 00245 bool computed_sparsity_already = 00246 ((_n_nz && !_n_nz->empty()) || 00247 (_n_oz && !_n_oz->empty())); 00248 this->comm().max(computed_sparsity_already); 00249 if (computed_sparsity_already && 00250 matrix.need_full_sparsity_pattern()) 00251 { 00252 // We'd better have already computed the full sparsity pattern 00253 // if we need it here 00254 libmesh_assert(need_full_sparsity_pattern); 00255 libmesh_assert(_sp.get()); 00256 00257 matrix.update_sparsity_pattern (_sp->sparsity_pattern); 00258 } 00259 00260 if (matrix.need_full_sparsity_pattern()) 00261 need_full_sparsity_pattern = true; 00262 } 00263 00264 00265 00266 bool DofMap::is_attached (SparseMatrix<Number>& matrix) 00267 { 00268 return (std::find(_matrices.begin(), _matrices.end(), 00269 &matrix) != _matrices.end()); 00270 } 00271 00272 00273 00274 DofObject* DofMap::node_ptr(MeshBase& mesh, dof_id_type i) const 00275 { 00276 return mesh.node_ptr(i); 00277 } 00278 00279 00280 00281 DofObject* DofMap::elem_ptr(MeshBase& mesh, dof_id_type i) const 00282 { 00283 return mesh.elem(i); 00284 } 00285 00286 00287 00288 template <typename iterator_type> 00289 void DofMap::set_nonlocal_dof_objects(iterator_type objects_begin, 00290 iterator_type objects_end, 00291 MeshBase &mesh, 00292 dofobject_accessor objects) 00293 { 00294 // This function must be run on all processors at once 00295 parallel_object_only(); 00296 00297 // First, iterate over local objects to find out how many 00298 // are on each processor 00299 std::vector<dof_id_type> 00300 ghost_objects_from_proc(this->n_processors(), 0); 00301 00302 iterator_type it = objects_begin; 00303 00304 for (; it != objects_end; ++it) 00305 { 00306 DofObject *obj = *it; 00307 00308 if (obj) 00309 { 00310 processor_id_type obj_procid = obj->processor_id(); 00311 // We'd better be completely partitioned by now 00312 libmesh_assert_not_equal_to (obj_procid, DofObject::invalid_processor_id); 00313 ghost_objects_from_proc[obj_procid]++; 00314 } 00315 } 00316 00317 std::vector<dof_id_type> objects_on_proc(this->n_processors(), 0); 00318 this->comm().allgather(ghost_objects_from_proc[this->processor_id()], 00319 objects_on_proc); 00320 00321 #ifdef DEBUG 00322 for (processor_id_type p=0; p != this->n_processors(); ++p) 00323 libmesh_assert_less_equal (ghost_objects_from_proc[p], objects_on_proc[p]); 00324 #endif 00325 00326 // Request sets to send to each processor 00327 std::vector<std::vector<dof_id_type> > 00328 requested_ids(this->n_processors()); 00329 00330 // We know how many of our objects live on each processor, so 00331 // reserve() space for requests from each. 00332 for (processor_id_type p=0; p != this->n_processors(); ++p) 00333 if (p != this->processor_id()) 00334 requested_ids[p].reserve(ghost_objects_from_proc[p]); 00335 00336 for (it = objects_begin; it != objects_end; ++it) 00337 { 00338 DofObject *obj = *it; 00339 if (obj->processor_id() != DofObject::invalid_processor_id) 00340 requested_ids[obj->processor_id()].push_back(obj->id()); 00341 } 00342 #ifdef DEBUG 00343 for (processor_id_type p=0; p != this->n_processors(); ++p) 00344 libmesh_assert_equal_to (requested_ids[p].size(), ghost_objects_from_proc[p]); 00345 #endif 00346 00347 // Next set ghost object n_comps from other processors 00348 for (processor_id_type p=1; p != this->n_processors(); ++p) 00349 { 00350 // Trade my requests with processor procup and procdown 00351 processor_id_type procup = 00352 cast_int<processor_id_type>((this->processor_id() + p) % 00353 this->n_processors()); 00354 processor_id_type procdown = 00355 cast_int<processor_id_type>((this->n_processors() + 00356 this->processor_id() - p) % 00357 this->n_processors()); 00358 std::vector<dof_id_type> request_to_fill; 00359 this->comm().send_receive(procup, requested_ids[procup], 00360 procdown, request_to_fill); 00361 00362 // Fill those requests 00363 const unsigned int 00364 sys_num = this->sys_number(), 00365 n_var_groups = this->n_variable_groups(); 00366 00367 std::vector<dof_id_type> ghost_data 00368 (request_to_fill.size() * 2 * n_var_groups); 00369 00370 for (std::size_t i=0; i != request_to_fill.size(); ++i) 00371 { 00372 DofObject *requested = (this->*objects)(mesh, request_to_fill[i]); 00373 libmesh_assert(requested); 00374 libmesh_assert_equal_to (requested->processor_id(), this->processor_id()); 00375 libmesh_assert_equal_to (requested->n_var_groups(sys_num), n_var_groups); 00376 for (unsigned int vg=0; vg != n_var_groups; ++vg) 00377 { 00378 unsigned int n_comp_g = 00379 requested->n_comp_group(sys_num, vg); 00380 ghost_data[i*2*n_var_groups+vg] = n_comp_g; 00381 dof_id_type my_first_dof = n_comp_g ? 00382 requested->vg_dof_base(sys_num, vg) : 0; 00383 libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id); 00384 ghost_data[i*2*n_var_groups+n_var_groups+vg] = my_first_dof; 00385 } 00386 } 00387 00388 // Trade back the results 00389 std::vector<dof_id_type> filled_request; 00390 this->comm().send_receive(procdown, ghost_data, 00391 procup, filled_request); 00392 00393 // And copy the id changes we've now been informed of 00394 libmesh_assert_equal_to (filled_request.size(), 00395 requested_ids[procup].size() * 2 * n_var_groups); 00396 for (std::size_t i=0; i != requested_ids[procup].size(); ++i) 00397 { 00398 DofObject *requested = (this->*objects)(mesh, requested_ids[procup][i]); 00399 libmesh_assert(requested); 00400 libmesh_assert_equal_to (requested->processor_id(), procup); 00401 for (unsigned int vg=0; vg != n_var_groups; ++vg) 00402 { 00403 unsigned int n_comp_g = 00404 cast_int<unsigned int>(filled_request[i*2*n_var_groups+vg]); 00405 requested->set_n_comp_group(sys_num, vg, n_comp_g); 00406 if (n_comp_g) 00407 { 00408 dof_id_type my_first_dof = 00409 filled_request[i*2*n_var_groups+n_var_groups+vg]; 00410 libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id); 00411 requested->set_vg_dof_base 00412 (sys_num, vg, my_first_dof); 00413 } 00414 } 00415 } 00416 } 00417 00418 #ifdef DEBUG 00419 // Double check for invalid dofs 00420 for (it = objects_begin; it != objects_end; ++it) 00421 { 00422 DofObject *obj = *it; 00423 libmesh_assert (obj); 00424 unsigned int num_variables = obj->n_vars(this->sys_number()); 00425 for (unsigned int v=0; v != num_variables; ++v) 00426 { 00427 unsigned int n_comp = 00428 obj->n_comp(this->sys_number(), v); 00429 dof_id_type my_first_dof = n_comp ? 00430 obj->dof_number(this->sys_number(), v, 0) : 0; 00431 libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id); 00432 } 00433 } 00434 #endif 00435 } 00436 00437 00438 00439 void DofMap::reinit(MeshBase& mesh) 00440 { 00441 libmesh_assert (mesh.is_prepared()); 00442 00443 START_LOG("reinit()", "DofMap"); 00444 00445 const unsigned int 00446 sys_num = this->sys_number(), 00447 n_var_groups = this->n_variable_groups(); 00448 00449 // The DofObjects need to know how many variable groups we have, and 00450 // how many variables there are in each group. 00451 std::vector<unsigned int> n_vars_per_group; n_vars_per_group.reserve (n_var_groups); 00452 00453 for (unsigned int vg=0; vg<n_var_groups; vg++) 00454 n_vars_per_group.push_back (this->variable_group(vg).n_variables()); 00455 00456 #ifdef LIBMESH_ENABLE_AMR 00457 00458 //------------------------------------------------------------ 00459 // Clear the old_dof_objects for all the nodes 00460 // and elements so that we can overwrite them 00461 { 00462 MeshBase::node_iterator node_it = mesh.nodes_begin(); 00463 const MeshBase::node_iterator node_end = mesh.nodes_end(); 00464 00465 for ( ; node_it != node_end; ++node_it) 00466 { 00467 (*node_it)->clear_old_dof_object(); 00468 libmesh_assert (!(*node_it)->old_dof_object); 00469 } 00470 00471 MeshBase::element_iterator elem_it = mesh.elements_begin(); 00472 const MeshBase::element_iterator elem_end = mesh.elements_end(); 00473 00474 for ( ; elem_it != elem_end; ++elem_it) 00475 { 00476 (*elem_it)->clear_old_dof_object(); 00477 libmesh_assert (!(*elem_it)->old_dof_object); 00478 } 00479 } 00480 00481 00482 //------------------------------------------------------------ 00483 // Set the old_dof_objects for the elements that 00484 // weren't just created, if these old dof objects 00485 // had variables 00486 { 00487 MeshBase::element_iterator elem_it = mesh.elements_begin(); 00488 const MeshBase::element_iterator elem_end = mesh.elements_end(); 00489 00490 for ( ; elem_it != elem_end; ++elem_it) 00491 { 00492 Elem* elem = *elem_it; 00493 00494 // Skip the elements that were just refined 00495 if (elem->refinement_flag() == Elem::JUST_REFINED) continue; 00496 00497 for (unsigned int n=0; n<elem->n_nodes(); n++) 00498 { 00499 Node* node = elem->get_node(n); 00500 00501 if (node->old_dof_object == NULL) 00502 if (node->has_dofs(sys_num)) 00503 node->set_old_dof_object(); 00504 } 00505 00506 libmesh_assert (!elem->old_dof_object); 00507 00508 if (elem->has_dofs(sys_num)) 00509 elem->set_old_dof_object(); 00510 } 00511 } 00512 00513 #endif // #ifdef LIBMESH_ENABLE_AMR 00514 00515 00516 //------------------------------------------------------------ 00517 // Then set the number of variables for each \p DofObject 00518 // equal to n_variables() for this system. This will 00519 // handle new \p DofObjects that may have just been created 00520 { 00521 // All the nodes 00522 MeshBase::node_iterator node_it = mesh.nodes_begin(); 00523 const MeshBase::node_iterator node_end = mesh.nodes_end(); 00524 00525 for ( ; node_it != node_end; ++node_it) 00526 (*node_it)->set_n_vars_per_group(sys_num, n_vars_per_group); 00527 00528 // All the elements 00529 MeshBase::element_iterator elem_it = mesh.elements_begin(); 00530 const MeshBase::element_iterator elem_end = mesh.elements_end(); 00531 00532 for ( ; elem_it != elem_end; ++elem_it) 00533 (*elem_it)->set_n_vars_per_group(sys_num, n_vars_per_group); 00534 } 00535 00536 00537 // Zero _n_SCALAR_dofs, it will be updated below. 00538 this->_n_SCALAR_dofs = 0; 00539 00540 //------------------------------------------------------------ 00541 // Next allocate space for the DOF indices 00542 for (unsigned int vg=0; vg<n_var_groups; vg++) 00543 { 00544 const VariableGroup &vg_description = this->variable_group(vg); 00545 00546 const unsigned int n_var_in_group = vg_description.n_variables(); 00547 const FEType& base_fe_type = vg_description.type(); 00548 00549 // Don't need to loop over elements for a SCALAR variable 00550 // Just increment _n_SCALAR_dofs 00551 if(base_fe_type.family == SCALAR) 00552 { 00553 this->_n_SCALAR_dofs += base_fe_type.order*n_var_in_group; 00554 continue; 00555 } 00556 00557 // This should be constant even on p-refined elements 00558 const bool extra_hanging_dofs = 00559 FEInterface::extra_hanging_dofs(base_fe_type); 00560 00561 // For all the active elements 00562 MeshBase::element_iterator elem_it = mesh.active_elements_begin(); 00563 const MeshBase::element_iterator elem_end = mesh.active_elements_end(); 00564 00565 // Count vertex degrees of freedom first 00566 for ( ; elem_it != elem_end; ++elem_it) 00567 { 00568 Elem* elem = *elem_it; 00569 libmesh_assert(elem); 00570 00571 // Skip the numbering if this variable is 00572 // not active on this element's subdomain 00573 if (!vg_description.active_on_subdomain(elem->subdomain_id())) 00574 continue; 00575 00576 const ElemType type = elem->type(); 00577 const unsigned int dim = elem->dim(); 00578 00579 FEType fe_type = base_fe_type; 00580 00581 #ifdef LIBMESH_ENABLE_AMR 00582 // Make sure we haven't done more p refinement than we can 00583 // handle 00584 if (elem->p_level() + base_fe_type.order > 00585 FEInterface::max_order(base_fe_type, type)) 00586 { 00587 # ifdef DEBUG 00588 if (FEInterface::max_order(base_fe_type,type) < static_cast<unsigned int>(base_fe_type.order)) 00589 libmesh_error_msg("ERROR: Finite element " \ 00590 << Utility::enum_to_string(base_fe_type.family) \ 00591 << " on geometric element " \ 00592 << Utility::enum_to_string(type) \ 00593 << "\nonly supports FEInterface::max_order = " \ 00594 << FEInterface::max_order(base_fe_type,type) \ 00595 << ", not fe_type.order = " \ 00596 << base_fe_type.order); 00597 00598 libMesh::err 00599 << "WARNING: Finite element " 00600 << Utility::enum_to_string(base_fe_type.family) 00601 << " on geometric element " 00602 << Utility::enum_to_string(type) << std::endl 00603 << "could not be p refined past FEInterface::max_order = " 00604 << FEInterface::max_order(base_fe_type,type) 00605 << std::endl; 00606 # endif 00607 elem->set_p_level(FEInterface::max_order(base_fe_type,type) 00608 - base_fe_type.order); 00609 } 00610 #endif 00611 00612 fe_type.order = static_cast<Order>(fe_type.order + 00613 elem->p_level()); 00614 00615 // Allocate the vertex DOFs 00616 for (unsigned int n=0; n<elem->n_nodes(); n++) 00617 { 00618 Node* node = elem->get_node(n); 00619 00620 if (elem->is_vertex(n)) 00621 { 00622 const unsigned int old_node_dofs = 00623 node->n_comp_group(sys_num, vg); 00624 00625 const unsigned int vertex_dofs = 00626 std::max(FEInterface::n_dofs_at_node(dim, fe_type, 00627 type, n), 00628 old_node_dofs); 00629 00630 // Some discontinuous FEs have no vertex dofs 00631 if (vertex_dofs > old_node_dofs) 00632 { 00633 node->set_n_comp_group(sys_num, vg, 00634 vertex_dofs); 00635 00636 // Abusing dof_number to set a "this is a 00637 // vertex" flag 00638 node->set_vg_dof_base(sys_num, vg, 00639 vertex_dofs); 00640 00641 // libMesh::out << "sys_num,vg,old_node_dofs,vertex_dofs=" 00642 // << sys_num << "," 00643 // << vg << "," 00644 // << old_node_dofs << "," 00645 // << vertex_dofs << '\n', 00646 // node->debug_buffer(); 00647 00648 // libmesh_assert_equal_to (vertex_dofs, node->n_comp(sys_num, vg)); 00649 // libmesh_assert_equal_to (vertex_dofs, node->vg_dof_base(sys_num, vg)); 00650 } 00651 } 00652 } 00653 } // done counting vertex dofs 00654 00655 // count edge & face dofs next 00656 elem_it = mesh.active_elements_begin(); 00657 00658 for ( ; elem_it != elem_end; ++elem_it) 00659 { 00660 Elem* elem = *elem_it; 00661 libmesh_assert(elem); 00662 00663 // Skip the numbering if this variable is 00664 // not active on this element's subdomain 00665 if (!vg_description.active_on_subdomain(elem->subdomain_id())) 00666 continue; 00667 00668 const ElemType type = elem->type(); 00669 const unsigned int dim = elem->dim(); 00670 00671 FEType fe_type = base_fe_type; 00672 fe_type.order = static_cast<Order>(fe_type.order + 00673 elem->p_level()); 00674 00675 // Allocate the edge and face DOFs 00676 for (unsigned int n=0; n<elem->n_nodes(); n++) 00677 { 00678 Node* node = elem->get_node(n); 00679 00680 const unsigned int old_node_dofs = 00681 node->n_comp_group(sys_num, vg); 00682 00683 const unsigned int vertex_dofs = old_node_dofs? 00684 cast_int<unsigned int>(node->vg_dof_base (sys_num,vg)):0; 00685 00686 const unsigned int new_node_dofs = 00687 FEInterface::n_dofs_at_node(dim, fe_type, type, n); 00688 00689 // We've already allocated vertex DOFs 00690 if (elem->is_vertex(n)) 00691 { 00692 libmesh_assert_greater_equal (old_node_dofs, vertex_dofs); 00693 // //if (vertex_dofs < new_node_dofs) 00694 // libMesh::out << "sys_num,vg,old_node_dofs,vertex_dofs,new_node_dofs=" 00695 // << sys_num << "," 00696 // << vg << "," 00697 // << old_node_dofs << "," 00698 // << vertex_dofs << "," 00699 // << new_node_dofs << '\n', 00700 // node->debug_buffer(); 00701 00702 libmesh_assert_greater_equal (vertex_dofs, new_node_dofs); 00703 } 00704 // We need to allocate the rest 00705 else 00706 { 00707 // If this has no dofs yet, it needs no vertex 00708 // dofs, so we just give it edge or face dofs 00709 if (!old_node_dofs) 00710 { 00711 node->set_n_comp_group(sys_num, vg, 00712 new_node_dofs); 00713 // Abusing dof_number to set a "this has no 00714 // vertex dofs" flag 00715 if (new_node_dofs) 00716 node->set_vg_dof_base(sys_num, vg, 00717 0); 00718 } 00719 00720 // If this has dofs, but has no vertex dofs, 00721 // it may still need more edge or face dofs if 00722 // we're p-refined. 00723 else if (vertex_dofs == 0) 00724 { 00725 if (new_node_dofs > old_node_dofs) 00726 { 00727 node->set_n_comp_group(sys_num, vg, 00728 new_node_dofs); 00729 00730 node->set_vg_dof_base(sys_num, vg, 00731 vertex_dofs); 00732 } 00733 } 00734 // If this is another element's vertex, 00735 // add more (non-overlapping) edge/face dofs if 00736 // necessary 00737 else if (extra_hanging_dofs) 00738 { 00739 if (new_node_dofs > old_node_dofs - vertex_dofs) 00740 { 00741 node->set_n_comp_group(sys_num, vg, 00742 vertex_dofs + new_node_dofs); 00743 00744 node->set_vg_dof_base(sys_num, vg, 00745 vertex_dofs); 00746 } 00747 } 00748 // If this is another element's vertex, add any 00749 // (overlapping) edge/face dofs if necessary 00750 else 00751 { 00752 libmesh_assert_greater_equal (old_node_dofs, vertex_dofs); 00753 if (new_node_dofs > old_node_dofs) 00754 { 00755 node->set_n_comp_group(sys_num, vg, 00756 new_node_dofs); 00757 00758 node->set_vg_dof_base (sys_num, vg, 00759 vertex_dofs); 00760 } 00761 } 00762 } 00763 } 00764 // Allocate the element DOFs 00765 const unsigned int dofs_per_elem = 00766 FEInterface::n_dofs_per_elem(dim, fe_type, 00767 type); 00768 00769 elem->set_n_comp_group(sys_num, vg, dofs_per_elem); 00770 00771 } 00772 } // end loop over variable groups 00773 00774 // Calling DofMap::reinit() by itself makes little sense, 00775 // so we won't bother with nonlocal DofObjects. 00776 // Those will be fixed by distribute_dofs 00777 00778 //------------------------------------------------------------ 00779 // Finally, clear all the current DOF indices 00780 // (distribute_dofs expects them cleared!) 00781 this->invalidate_dofs(mesh); 00782 00783 STOP_LOG("reinit()", "DofMap"); 00784 } 00785 00786 00787 00788 void DofMap::invalidate_dofs(MeshBase& mesh) const 00789 { 00790 const unsigned int sys_num = this->sys_number(); 00791 00792 // All the nodes 00793 MeshBase::node_iterator node_it = mesh.nodes_begin(); 00794 const MeshBase::node_iterator node_end = mesh.nodes_end(); 00795 00796 for ( ; node_it != node_end; ++node_it) 00797 (*node_it)->invalidate_dofs(sys_num); 00798 00799 // All the elements 00800 MeshBase::element_iterator elem_it = mesh.active_elements_begin(); 00801 const MeshBase::element_iterator elem_end = mesh.active_elements_end(); 00802 00803 for ( ; elem_it != elem_end; ++elem_it) 00804 (*elem_it)->invalidate_dofs(sys_num); 00805 } 00806 00807 00808 00809 void DofMap::clear() 00810 { 00811 // we don't want to clear 00812 // the coupling matrix! 00813 // It should not change... 00814 //_dof_coupling->clear(); 00815 00816 _variables.clear(); 00817 _variable_groups.clear(); 00818 _first_df.clear(); 00819 _end_df.clear(); 00820 _first_scalar_df.clear(); 00821 _send_list.clear(); 00822 this->clear_sparsity(); 00823 need_full_sparsity_pattern = false; 00824 00825 #ifdef LIBMESH_ENABLE_AMR 00826 00827 _dof_constraints.clear(); 00828 _stashed_dof_constraints.clear(); 00829 _primal_constraint_values.clear(); 00830 _adjoint_constraint_values.clear(); 00831 _n_old_dfs = 0; 00832 _first_old_df.clear(); 00833 _end_old_df.clear(); 00834 _first_old_scalar_df.clear(); 00835 00836 #endif 00837 00838 _matrices.clear(); 00839 00840 _n_dfs = 0; 00841 } 00842 00843 00844 00845 void DofMap::distribute_dofs (MeshBase& mesh) 00846 { 00847 // This function must be run on all processors at once 00848 parallel_object_only(); 00849 00850 // Log how long it takes to distribute the degrees of freedom 00851 START_LOG("distribute_dofs()", "DofMap"); 00852 00853 libmesh_assert (mesh.is_prepared()); 00854 00855 const processor_id_type proc_id = this->processor_id(); 00856 const processor_id_type n_proc = this->n_processors(); 00857 00858 // libmesh_assert_greater (this->n_variables(), 0); 00859 libmesh_assert_less (proc_id, n_proc); 00860 00861 // re-init in case the mesh has changed 00862 this->reinit(mesh); 00863 00864 // By default distribute variables in a 00865 // var-major fashion, but allow run-time 00866 // specification 00867 bool node_major_dofs = libMesh::on_command_line ("--node_major_dofs"); 00868 00869 // The DOF counter, will be incremented as we encounter 00870 // new degrees of freedom 00871 dof_id_type next_free_dof = 0; 00872 00873 // Clear the send list before we rebuild it 00874 _send_list.clear(); 00875 00876 // Set temporary DOF indices on this processor 00877 if (node_major_dofs) 00878 this->distribute_local_dofs_node_major (next_free_dof, mesh); 00879 else 00880 this->distribute_local_dofs_var_major (next_free_dof, mesh); 00881 00882 // Get DOF counts on all processors 00883 std::vector<dof_id_type> dofs_on_proc(n_proc, 0); 00884 this->comm().allgather(next_free_dof, dofs_on_proc); 00885 00886 // Resize and fill the _first_df and _end_df arrays 00887 #ifdef LIBMESH_ENABLE_AMR 00888 _first_old_df = _first_df; 00889 _end_old_df = _end_df; 00890 #endif 00891 00892 _first_df.resize(n_proc); 00893 _end_df.resize (n_proc); 00894 00895 // Get DOF offsets 00896 _first_df[0] = 0; 00897 for (processor_id_type i=1; i < n_proc; ++i) 00898 _first_df[i] = _end_df[i-1] = _first_df[i-1] + dofs_on_proc[i-1]; 00899 _end_df[n_proc-1] = _first_df[n_proc-1] + dofs_on_proc[n_proc-1]; 00900 00901 // Clear all the current DOF indices 00902 // (distribute_dofs expects them cleared!) 00903 this->invalidate_dofs(mesh); 00904 00905 next_free_dof = _first_df[proc_id]; 00906 00907 // Set permanent DOF indices on this processor 00908 if (node_major_dofs) 00909 this->distribute_local_dofs_node_major (next_free_dof, mesh); 00910 else 00911 this->distribute_local_dofs_var_major (next_free_dof, mesh); 00912 00913 libmesh_assert_equal_to (next_free_dof, _end_df[proc_id]); 00914 00915 //------------------------------------------------------------ 00916 // At this point, all n_comp and dof_number values on local 00917 // DofObjects should be correct, but a ParallelMesh might have 00918 // incorrect values on non-local DofObjects. Let's request the 00919 // correct values from each other processor. 00920 00921 if (this->n_processors() > 1) 00922 { 00923 this->set_nonlocal_dof_objects(mesh.nodes_begin(), 00924 mesh.nodes_end(), 00925 mesh, &DofMap::node_ptr); 00926 00927 this->set_nonlocal_dof_objects(mesh.elements_begin(), 00928 mesh.elements_end(), 00929 mesh, &DofMap::elem_ptr); 00930 } 00931 00932 #ifdef DEBUG 00933 { 00934 const unsigned int 00935 sys_num = this->sys_number(); 00936 00937 // Processors should all agree on DoF ids 00938 MeshTools::libmesh_assert_valid_dof_ids(mesh); 00939 00940 // DoF processor ids should match DofObject processor ids 00941 MeshBase::const_node_iterator node_it = mesh.nodes_begin(); 00942 const MeshBase::const_node_iterator node_end = mesh.nodes_end(); 00943 for ( ; node_it != node_end; ++node_it) 00944 { 00945 DofObject const * const dofobj = *node_it; 00946 const processor_id_type obj_proc_id = dofobj->processor_id(); 00947 00948 for (unsigned int v=0; v != dofobj->n_vars(sys_num); ++v) 00949 for (unsigned int c=0; c != dofobj->n_comp(sys_num,v); ++c) 00950 { 00951 const dof_id_type dofid = dofobj->dof_number(sys_num,v,c); 00952 libmesh_assert_greater_equal (dofid, this->first_dof(obj_proc_id)); 00953 libmesh_assert_less (dofid, this->end_dof(obj_proc_id)); 00954 } 00955 } 00956 00957 MeshBase::const_element_iterator elem_it = mesh.elements_begin(); 00958 const MeshBase::const_element_iterator elem_end = mesh.elements_end(); 00959 for ( ; elem_it != elem_end; ++elem_it) 00960 { 00961 DofObject const * const dofobj = *elem_it; 00962 const processor_id_type obj_proc_id = dofobj->processor_id(); 00963 00964 for (unsigned int v=0; v != dofobj->n_vars(sys_num); ++v) 00965 for (unsigned int c=0; c != dofobj->n_comp(sys_num,v); ++c) 00966 { 00967 const dof_id_type dofid = dofobj->dof_number(sys_num,v,c); 00968 libmesh_assert_greater_equal (dofid, this->first_dof(obj_proc_id)); 00969 libmesh_assert_less (dofid, this->end_dof(obj_proc_id)); 00970 } 00971 } 00972 } 00973 #endif 00974 00975 // Set the total number of degrees of freedom, then start finding 00976 // SCALAR degrees of freedom 00977 #ifdef LIBMESH_ENABLE_AMR 00978 _n_old_dfs = _n_dfs; 00979 _first_old_scalar_df = _first_scalar_df; 00980 #endif 00981 _n_dfs = _end_df[n_proc-1]; 00982 _first_scalar_df.clear(); 00983 _first_scalar_df.resize(this->n_variables(), DofObject::invalid_id); 00984 dof_id_type current_SCALAR_dof_index = n_dofs() - n_SCALAR_dofs(); 00985 00986 // Calculate and cache the initial DoF indices for SCALAR variables. 00987 // This is an O(N_vars) calculation so we want to do it once per 00988 // renumbering rather than once per SCALAR_dof_indices() call 00989 00990 for (unsigned int v=0; v<this->n_variables(); v++) 00991 if(this->variable(v).type().family == SCALAR) 00992 { 00993 _first_scalar_df[v] = current_SCALAR_dof_index; 00994 current_SCALAR_dof_index += this->variable(v).type().order; 00995 } 00996 00997 STOP_LOG("distribute_dofs()", "DofMap"); 00998 00999 // Note that in the add_neighbors_to_send_list nodes on processor 01000 // boundaries that are shared by multiple elements are added for 01001 // each element. 01002 this->add_neighbors_to_send_list(mesh); 01003 01004 // Here we used to clean up that data structure; now System and 01005 // EquationSystems call that for us, after we've added constraint 01006 // dependencies to the send_list too. 01007 // this->sort_send_list (); 01008 } 01009 01010 01011 void DofMap::local_variable_indices(std::vector<dof_id_type>& idx, 01012 const MeshBase& mesh, 01013 unsigned int var_num) const 01014 { 01015 const unsigned int sys_num = this->sys_number(); 01016 01017 // If this isn't a SCALAR variable, we need to find all its field 01018 // dofs on the mesh 01019 if (this->variable_type(var_num).family != SCALAR) 01020 { 01021 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 01022 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 01023 01024 for ( ; elem_it != elem_end; ++elem_it) 01025 { 01026 // Only count dofs connected to active 01027 // elements on this processor. 01028 Elem* elem = *elem_it; 01029 const unsigned int n_nodes = elem->n_nodes(); 01030 01031 // First get any new nodal DOFS 01032 for (unsigned int n=0; n<n_nodes; n++) 01033 { 01034 Node* node = elem->get_node(n); 01035 01036 if (node->processor_id() < this->processor_id()) 01037 continue; 01038 01039 const unsigned int n_comp = node->n_comp(sys_num, var_num); 01040 for(unsigned int i=0; i<n_comp; i++) 01041 { 01042 const dof_id_type index = node->dof_number(sys_num,var_num,i); 01043 libmesh_assert_greater_equal (index, this->first_dof()); 01044 libmesh_assert_less (index, this->end_dof()); 01045 01046 if (idx.empty() || index > idx.back()) 01047 idx.push_back(index); 01048 } 01049 } 01050 01051 // Next get any new element DOFS 01052 const unsigned int n_comp = elem->n_comp(sys_num, var_num); 01053 for(unsigned int i=0; i<n_comp; i++) 01054 { 01055 const dof_id_type index = elem->dof_number(sys_num,var_num,i); 01056 if (idx.empty() || index > idx.back()) 01057 idx.push_back(index); 01058 } 01059 } // done looping over elements 01060 01061 01062 // we may have missed assigning DOFs to nodes that we own 01063 // but to which we have no connected elements matching our 01064 // variable restriction criterion. this will happen, for example, 01065 // if variable V is restricted to subdomain S. We may not own 01066 // any elements which live in S, but we may own nodes which are 01067 // *connected* to elements which do. in this scenario these nodes 01068 // will presently have unnumbered DOFs. we need to take care of 01069 // them here since we own them and no other processor will touch them. 01070 { 01071 MeshBase::const_node_iterator node_it = mesh.local_nodes_begin(); 01072 const MeshBase::const_node_iterator node_end = mesh.local_nodes_end(); 01073 01074 for (; node_it != node_end; ++node_it) 01075 { 01076 Node *node = *node_it; 01077 libmesh_assert(node); 01078 01079 const unsigned int n_comp = node->n_comp(sys_num, var_num); 01080 for(unsigned int i=0; i<n_comp; i++) 01081 { 01082 const dof_id_type index = node->dof_number(sys_num,var_num,i); 01083 if (idx.empty() || index > idx.back()) 01084 idx.push_back(index); 01085 } 01086 } 01087 } 01088 } 01089 // Otherwise, count up the SCALAR dofs, if we're on the processor 01090 // that holds this SCALAR variable 01091 else if ( this->processor_id() == (this->n_processors()-1) ) 01092 { 01093 std::vector<dof_id_type> di_scalar; 01094 this->SCALAR_dof_indices(di_scalar,var_num); 01095 idx.insert( idx.end(), di_scalar.begin(), di_scalar.end()); 01096 } 01097 } 01098 01099 01100 void DofMap::distribute_local_dofs_node_major(dof_id_type &next_free_dof, 01101 MeshBase& mesh) 01102 { 01103 const unsigned int sys_num = this->sys_number(); 01104 const unsigned int n_var_groups = this->n_variable_groups(); 01105 01106 //------------------------------------------------------------------------- 01107 // First count and assign temporary numbers to local dofs 01108 MeshBase::element_iterator elem_it = mesh.active_local_elements_begin(); 01109 const MeshBase::element_iterator elem_end = mesh.active_local_elements_end(); 01110 01111 for ( ; elem_it != elem_end; ++elem_it) 01112 { 01113 // Only number dofs connected to active 01114 // elements on this processor. 01115 Elem* elem = *elem_it; 01116 const unsigned int n_nodes = elem->n_nodes(); 01117 01118 // First number the nodal DOFS 01119 for (unsigned int n=0; n<n_nodes; n++) 01120 { 01121 Node* node = elem->get_node(n); 01122 01123 for (unsigned vg=0; vg<n_var_groups; vg++) 01124 { 01125 const VariableGroup &vg_description(this->variable_group(vg)); 01126 01127 if( (vg_description.type().family != SCALAR) && 01128 (vg_description.active_on_subdomain(elem->subdomain_id())) ) 01129 { 01130 // assign dof numbers (all at once) if this is 01131 // our node and if they aren't already there 01132 if ((node->n_comp_group(sys_num,vg) > 0) && 01133 (node->processor_id() == this->processor_id()) && 01134 (node->vg_dof_base(sys_num,vg) == 01135 DofObject::invalid_id)) 01136 { 01137 node->set_vg_dof_base(sys_num, 01138 vg, 01139 next_free_dof); 01140 next_free_dof += (vg_description.n_variables()* 01141 node->n_comp_group(sys_num,vg)); 01142 //node->debug_buffer(); 01143 } 01144 } 01145 } 01146 } 01147 01148 // Now number the element DOFS 01149 for (unsigned vg=0; vg<n_var_groups; vg++) 01150 { 01151 const VariableGroup &vg_description(this->variable_group(vg)); 01152 01153 if ( (vg_description.type().family != SCALAR) && 01154 (vg_description.active_on_subdomain(elem->subdomain_id())) ) 01155 if (elem->n_comp_group(sys_num,vg) > 0) 01156 { 01157 libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg), 01158 DofObject::invalid_id); 01159 01160 elem->set_vg_dof_base(sys_num, 01161 vg, 01162 next_free_dof); 01163 01164 next_free_dof += (vg_description.n_variables()* 01165 elem->n_comp(sys_num,vg)); 01166 } 01167 } 01168 } // done looping over elements 01169 01170 01171 // we may have missed assigning DOFs to nodes that we own 01172 // but to which we have no connected elements matching our 01173 // variable restriction criterion. this will happen, for example, 01174 // if variable V is restricted to subdomain S. We may not own 01175 // any elements which live in S, but we may own nodes which are 01176 // *connected* to elements which do. in this scenario these nodes 01177 // will presently have unnumbered DOFs. we need to take care of 01178 // them here since we own them and no other processor will touch them. 01179 { 01180 MeshBase::node_iterator node_it = mesh.local_nodes_begin(); 01181 const MeshBase::node_iterator node_end = mesh.local_nodes_end(); 01182 01183 for (; node_it != node_end; ++node_it) 01184 { 01185 Node *node = *node_it; 01186 libmesh_assert(node); 01187 01188 for (unsigned vg=0; vg<n_var_groups; vg++) 01189 { 01190 const VariableGroup &vg_description(this->variable_group(vg)); 01191 01192 if (node->n_comp_group(sys_num,vg)) 01193 if (node->vg_dof_base(sys_num,vg) == DofObject::invalid_id) 01194 { 01195 node->set_vg_dof_base (sys_num, 01196 vg, 01197 next_free_dof); 01198 01199 next_free_dof += (vg_description.n_variables()* 01200 node->n_comp(sys_num,vg)); 01201 } 01202 } 01203 } 01204 } 01205 01206 // Finally, count up the SCALAR dofs 01207 this->_n_SCALAR_dofs = 0; 01208 for (unsigned vg=0; vg<n_var_groups; vg++) 01209 { 01210 const VariableGroup &vg_description(this->variable_group(vg)); 01211 01212 if( vg_description.type().family == SCALAR ) 01213 { 01214 this->_n_SCALAR_dofs += (vg_description.n_variables()* 01215 vg_description.type().order); 01216 continue; 01217 } 01218 } 01219 01220 // Only increment next_free_dof if we're on the processor 01221 // that holds this SCALAR variable 01222 if ( this->processor_id() == (this->n_processors()-1) ) 01223 next_free_dof += _n_SCALAR_dofs; 01224 01225 #ifdef DEBUG 01226 { 01227 // libMesh::out << "next_free_dof=" << next_free_dof << std::endl 01228 // << "_n_SCALAR_dofs=" << _n_SCALAR_dofs << std::endl; 01229 01230 // Make sure we didn't miss any nodes 01231 MeshTools::libmesh_assert_valid_procids<Node>(mesh); 01232 01233 MeshBase::node_iterator node_it = mesh.local_nodes_begin(); 01234 const MeshBase::node_iterator node_end = mesh.local_nodes_end(); 01235 for (; node_it != node_end; ++node_it) 01236 { 01237 Node *obj = *node_it; 01238 libmesh_assert(obj); 01239 unsigned int n_var_g = obj->n_var_groups(this->sys_number()); 01240 for (unsigned int vg=0; vg != n_var_g; ++vg) 01241 { 01242 unsigned int n_comp_g = 01243 obj->n_comp_group(this->sys_number(), vg); 01244 dof_id_type my_first_dof = n_comp_g ? 01245 obj->vg_dof_base(this->sys_number(), vg) : 0; 01246 libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id); 01247 } 01248 } 01249 } 01250 #endif // DEBUG 01251 } 01252 01253 01254 01255 void DofMap::distribute_local_dofs_var_major(dof_id_type &next_free_dof, 01256 MeshBase& mesh) 01257 { 01258 const unsigned int sys_num = this->sys_number(); 01259 const unsigned int n_var_groups = this->n_variable_groups(); 01260 01261 //------------------------------------------------------------------------- 01262 // First count and assign temporary numbers to local dofs 01263 for (unsigned vg=0; vg<n_var_groups; vg++) 01264 { 01265 const VariableGroup &vg_description(this->variable_group(vg)); 01266 01267 const unsigned int n_vars_in_group = vg_description.n_variables(); 01268 01269 // Skip the SCALAR dofs 01270 if (vg_description.type().family == SCALAR) 01271 continue; 01272 01273 MeshBase::element_iterator elem_it = mesh.active_local_elements_begin(); 01274 const MeshBase::element_iterator elem_end = mesh.active_local_elements_end(); 01275 01276 for ( ; elem_it != elem_end; ++elem_it) 01277 { 01278 // Only number dofs connected to active 01279 // elements on this processor. 01280 Elem* elem = *elem_it; 01281 01282 // ... and only variables which are active on 01283 // on this element's subdomain 01284 if (!vg_description.active_on_subdomain(elem->subdomain_id())) 01285 continue; 01286 01287 const unsigned int n_nodes = elem->n_nodes(); 01288 01289 // First number the nodal DOFS 01290 for (unsigned int n=0; n<n_nodes; n++) 01291 { 01292 Node* node = elem->get_node(n); 01293 01294 // assign dof numbers (all at once) if this is 01295 // our node and if they aren't already there 01296 if ((node->n_comp_group(sys_num,vg) > 0) && 01297 (node->processor_id() == this->processor_id()) && 01298 (node->vg_dof_base(sys_num,vg) == 01299 DofObject::invalid_id)) 01300 { 01301 node->set_vg_dof_base(sys_num, 01302 vg, 01303 next_free_dof); 01304 01305 next_free_dof += (n_vars_in_group* 01306 node->n_comp_group(sys_num,vg)); 01307 } 01308 } 01309 01310 // Now number the element DOFS 01311 if (elem->n_comp_group(sys_num,vg) > 0) 01312 { 01313 libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg), 01314 DofObject::invalid_id); 01315 01316 elem->set_vg_dof_base(sys_num, 01317 vg, 01318 next_free_dof); 01319 01320 next_free_dof += (n_vars_in_group* 01321 elem->n_comp_group(sys_num,vg)); 01322 } 01323 } // end loop on elements 01324 01325 // we may have missed assigning DOFs to nodes that we own 01326 // but to which we have no connected elements matching our 01327 // variable restriction criterion. this will happen, for example, 01328 // if variable V is restricted to subdomain S. We may not own 01329 // any elements which live in S, but we may own nodes which are 01330 // *connected* to elements which do. in this scenario these nodes 01331 // will presently have unnumbered DOFs. we need to take care of 01332 // them here since we own them and no other processor will touch them. 01333 { 01334 MeshBase::node_iterator node_it = mesh.local_nodes_begin(); 01335 const MeshBase::node_iterator node_end = mesh.local_nodes_end(); 01336 01337 for (; node_it != node_end; ++node_it) 01338 { 01339 Node *node = *node_it; 01340 libmesh_assert(node); 01341 01342 if (node->n_comp_group(sys_num,vg)) 01343 if (node->vg_dof_base(sys_num,vg) == DofObject::invalid_id) 01344 { 01345 node->set_vg_dof_base (sys_num, 01346 vg, 01347 next_free_dof); 01348 01349 next_free_dof += (n_vars_in_group* 01350 node->n_comp_group(sys_num,vg)); 01351 } 01352 } 01353 } 01354 } // end loop on variable groups 01355 01356 // Finally, count up the SCALAR dofs 01357 this->_n_SCALAR_dofs = 0; 01358 for (unsigned vg=0; vg<n_var_groups; vg++) 01359 { 01360 const VariableGroup &vg_description(this->variable_group(vg)); 01361 01362 if( vg_description.type().family == SCALAR ) 01363 { 01364 this->_n_SCALAR_dofs += (vg_description.n_variables()* 01365 vg_description.type().order); 01366 continue; 01367 } 01368 } 01369 01370 // Only increment next_free_dof if we're on the processor 01371 // that holds this SCALAR variable 01372 if ( this->processor_id() == (this->n_processors()-1) ) 01373 next_free_dof += _n_SCALAR_dofs; 01374 01375 #ifdef DEBUG 01376 { 01377 // Make sure we didn't miss any nodes 01378 MeshTools::libmesh_assert_valid_procids<Node>(mesh); 01379 01380 MeshBase::node_iterator node_it = mesh.local_nodes_begin(); 01381 const MeshBase::node_iterator node_end = mesh.local_nodes_end(); 01382 for (; node_it != node_end; ++node_it) 01383 { 01384 Node *obj = *node_it; 01385 libmesh_assert(obj); 01386 unsigned int n_var_g = obj->n_var_groups(this->sys_number()); 01387 for (unsigned int vg=0; vg != n_var_g; ++vg) 01388 { 01389 unsigned int n_comp_g = 01390 obj->n_comp_group(this->sys_number(), vg); 01391 dof_id_type my_first_dof = n_comp_g ? 01392 obj->vg_dof_base(this->sys_number(), vg) : 0; 01393 libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id); 01394 } 01395 } 01396 } 01397 #endif // DEBUG 01398 } 01399 01400 01401 01402 void DofMap::add_neighbors_to_send_list(MeshBase& mesh) 01403 { 01404 START_LOG("add_neighbors_to_send_list()", "DofMap"); 01405 01406 const unsigned int sys_num = this->sys_number(); 01407 01408 //------------------------------------------------------------------------- 01409 // We need to add the DOFs from elements that live on neighboring processors 01410 // that are neighbors of the elements on the local processor 01411 //------------------------------------------------------------------------- 01412 01413 MeshBase::const_element_iterator local_elem_it 01414 = mesh.active_local_elements_begin(); 01415 const MeshBase::const_element_iterator local_elem_end 01416 = mesh.active_local_elements_end(); 01417 01418 std::vector<bool> node_on_processor(mesh.max_node_id(), false); 01419 std::vector<dof_id_type> di; 01420 std::vector<const Elem *> family; 01421 01422 // Loop over the active local elements, adding all active elements 01423 // that neighbor an active local element to the send list. 01424 for ( ; local_elem_it != local_elem_end; ++local_elem_it) 01425 { 01426 const Elem* elem = *local_elem_it; 01427 01428 // We may have non-local SCALAR dofs on a local element. 01429 // Normally we'd catch those on neighboring elements, but it's 01430 // possible that the neighbors are of different subdomains and 01431 // the SCALAR isn't supported on them. 01432 for (unsigned int v=0; v<this->n_variables(); v++) 01433 if(this->variable(v).type().family == SCALAR && 01434 this->variable(v).active_on_subdomain(elem->subdomain_id())) 01435 { 01436 // We asked for this variable, so add it to the vector. 01437 std::vector<dof_id_type> di_new; 01438 this->SCALAR_dof_indices(di_new,v); 01439 for (unsigned int i=0; i != di_new.size(); ++i) 01440 { 01441 const dof_id_type dof_index = di_new[i]; 01442 01443 if (dof_index < this->first_dof() || 01444 dof_index >= this->end_dof()) 01445 _send_list.push_back(dof_index); 01446 } 01447 } 01448 01449 // We may have non-local nodes on a local element. 01450 for (unsigned int n=0; n!=elem->n_nodes(); n++) 01451 { 01452 // Flag all the nodes of active local elements as seen, so 01453 // we can add nodal neighbor dofs to the send_list later. 01454 node_on_processor[elem->node(n)] = true; 01455 01456 // Add all remote dofs on these nodes to the send_list. 01457 // This is necessary in case those dofs are *not* also dofs 01458 // on neighbors; e.g. in the case of a HIERARCHIC's local 01459 // side which is only a vertex on the neighbor that owns it. 01460 const Node* node = elem->get_node(n); 01461 const unsigned n_vars = node->n_vars(sys_num); 01462 for (unsigned int v=0; v != n_vars; ++v) 01463 { 01464 const unsigned int n_comp = node->n_comp(sys_num, v); 01465 for (unsigned int c=0; c != n_comp; ++c) 01466 { 01467 const dof_id_type dof_index = node->dof_number(sys_num, v, c); 01468 if (dof_index < this->first_dof() || dof_index >= this->end_dof()) 01469 { 01470 _send_list.push_back(dof_index); 01471 // libmesh_here(); 01472 // libMesh::out << "sys_num,v,c,dof_index=" 01473 // << sys_num << ", " 01474 // << v << ", " 01475 // << c << ", " 01476 // << dof_index << '\n'; 01477 // node->debug_buffer(); 01478 } 01479 } 01480 } 01481 } 01482 01483 // Loop over the neighbors of those elements 01484 for (unsigned int s=0; s<elem->n_neighbors(); s++) 01485 if (elem->neighbor(s) != NULL) 01486 { 01487 family.clear(); 01488 01489 // Find all the active elements that neighbor elem 01490 #ifdef LIBMESH_ENABLE_AMR 01491 if (!elem->neighbor(s)->active()) 01492 elem->neighbor(s)->active_family_tree_by_neighbor(family, elem); 01493 else 01494 #endif 01495 family.push_back(elem->neighbor(s)); 01496 01497 for (dof_id_type i=0; i!=family.size(); ++i) 01498 // If the neighbor lives on a different processor 01499 if (family[i]->processor_id() != this->processor_id()) 01500 { 01501 // Get the DOF indices for this neighboring element 01502 this->dof_indices (family[i], di); 01503 01504 // Insert the remote DOF indices into the send list 01505 for (std::size_t j=0; j != di.size(); ++j) 01506 if (di[j] < this->first_dof() || 01507 di[j] >= this->end_dof()) 01508 _send_list.push_back(di[j]); 01509 } 01510 } 01511 } 01512 01513 // Now loop over all non_local active elements and add any missing 01514 // nodal-only neighbors. This will also get any dofs from nonlocal 01515 // nodes on local elements, because every nonlocal node exists on a 01516 // nonlocal nodal neighbor element. 01517 MeshBase::const_element_iterator elem_it 01518 = mesh.active_elements_begin(); 01519 const MeshBase::const_element_iterator elem_end 01520 = mesh.active_elements_end(); 01521 01522 for ( ; elem_it != elem_end; ++elem_it) 01523 { 01524 const Elem* elem = *elem_it; 01525 01526 // If this is one of our elements, we've already added it 01527 if (elem->processor_id() == this->processor_id()) 01528 continue; 01529 01530 // Do we need to add the element DOFs? 01531 bool add_elem_dofs = false; 01532 01533 // Check all the nodes of the element to see if it 01534 // shares a node with us 01535 for (unsigned int n=0; n!=elem->n_nodes(); n++) 01536 if (node_on_processor[elem->node(n)]) 01537 add_elem_dofs = true; 01538 01539 // Add the element degrees of freedom if it shares at 01540 // least one node. 01541 if (add_elem_dofs) 01542 { 01543 // Get the DOF indices for this neighboring element 01544 this->dof_indices (elem, di); 01545 01546 // Insert the remote DOF indices into the send list 01547 for (std::size_t j=0; j != di.size(); ++j) 01548 if (di[j] < this->first_dof() || 01549 di[j] >= this->end_dof()) 01550 _send_list.push_back(di[j]); 01551 } 01552 } 01553 01554 STOP_LOG("add_neighbors_to_send_list()", "DofMap"); 01555 } 01556 01557 01558 01559 void DofMap::prepare_send_list () 01560 { 01561 START_LOG("prepare_send_list()", "DofMap"); 01562 01563 // Check to see if we have any extra stuff to add to the send_list 01564 if (_extra_send_list_function) 01565 { 01566 if (_augment_send_list) 01567 { 01568 libmesh_here(); 01569 libMesh::out << "WARNING: You have specified both an extra send list function and object.\n" 01570 << " Are you sure this is what you meant to do??" 01571 << std::endl; 01572 } 01573 01574 _extra_send_list_function(_send_list, _extra_send_list_context); 01575 } 01576 01577 if (_augment_send_list) 01578 _augment_send_list->augment_send_list (_send_list); 01579 01580 // First sort the send list. After this 01581 // duplicated elements will be adjacent in the 01582 // vector 01583 std::sort(_send_list.begin(), _send_list.end()); 01584 01585 // Now use std::unique to remove duplicate entries 01586 std::vector<dof_id_type>::iterator new_end = 01587 std::unique (_send_list.begin(), _send_list.end()); 01588 01589 // Remove the end of the send_list. Use the "swap trick" 01590 // from Effective STL 01591 std::vector<dof_id_type> (_send_list.begin(), new_end).swap (_send_list); 01592 01593 STOP_LOG("prepare_send_list()", "DofMap"); 01594 } 01595 01596 01597 void DofMap::set_implicit_neighbor_dofs(bool implicit_neighbor_dofs) 01598 { 01599 _implicit_neighbor_dofs_initialized = true; 01600 _implicit_neighbor_dofs = implicit_neighbor_dofs; 01601 } 01602 01603 01604 bool DofMap::use_coupled_neighbor_dofs(const MeshBase& mesh) const 01605 { 01606 // If we were asked on the command line, then we need to 01607 // include sensitivities between neighbor degrees of freedom 01608 bool implicit_neighbor_dofs = 01609 libMesh::on_command_line ("--implicit_neighbor_dofs"); 01610 01611 // If the user specifies --implicit_neighbor_dofs 0, then 01612 // presumably he knows what he is doing and we won't try to 01613 // automatically turn it on even when all the variables are 01614 // discontinuous. 01615 if (implicit_neighbor_dofs) 01616 { 01617 // No flag provided defaults to 'true' 01618 int flag = 1; 01619 flag = libMesh::command_line_next ("--implicit_neighbor_dofs", flag); 01620 01621 if (!flag) 01622 { 01623 // The user said --implicit_neighbor_dofs 0, so he knows 01624 // what he is doing and really doesn't want it. 01625 return false; 01626 } 01627 } 01628 01629 // Possibly override the commandline option, if set_implicit_neighbor_dofs 01630 // has been called. 01631 if(_implicit_neighbor_dofs_initialized) 01632 { 01633 implicit_neighbor_dofs = _implicit_neighbor_dofs; 01634 01635 // Again, if the user explicitly says implicit_neighbor_dofs = false, 01636 // then we return here. 01637 if (!implicit_neighbor_dofs) 01638 return false; 01639 } 01640 01641 // Look at all the variables in this system. If every one is 01642 // discontinuous then the user must be doing DG/FVM, so be nice 01643 // and force implicit_neighbor_dofs=true. 01644 { 01645 bool all_discontinuous_dofs = true; 01646 01647 for (unsigned int var=0; var<this->n_variables(); var++) 01648 if (FEAbstract::build (mesh.mesh_dimension(), 01649 this->variable_type(var))->get_continuity() != DISCONTINUOUS) 01650 all_discontinuous_dofs = false; 01651 01652 if (all_discontinuous_dofs) 01653 implicit_neighbor_dofs = true; 01654 } 01655 01656 return implicit_neighbor_dofs; 01657 } 01658 01659 01660 01661 void DofMap::compute_sparsity(const MeshBase& mesh) 01662 { 01663 _sp = this->build_sparsity(mesh); 01664 01665 // It is possible that some \p SparseMatrix implementations want to 01666 // see it. Let them see it before we throw it away. 01667 std::vector<SparseMatrix<Number>* >::const_iterator 01668 pos = _matrices.begin(), 01669 end = _matrices.end(); 01670 01671 // If we need the full sparsity pattern, then we share a view of its 01672 // arrays, and we pass it in to the matrices. 01673 if (need_full_sparsity_pattern) 01674 { 01675 _n_nz = &_sp->n_nz; 01676 _n_oz = &_sp->n_oz; 01677 01678 for (; pos != end; ++pos) 01679 (*pos)->update_sparsity_pattern (_sp->sparsity_pattern); 01680 } 01681 // If we don't need the full sparsity pattern anymore, steal the 01682 // arrays we do need and free the rest of the memory 01683 else 01684 { 01685 if (!_n_nz) 01686 _n_nz = new std::vector<dof_id_type>(); 01687 _n_nz->swap(_sp->n_nz); 01688 if (!_n_oz) 01689 _n_oz = new std::vector<dof_id_type>(); 01690 _n_oz->swap(_sp->n_oz); 01691 01692 _sp.reset(); 01693 } 01694 } 01695 01696 01697 01698 void DofMap::clear_sparsity() 01699 { 01700 if (need_full_sparsity_pattern) 01701 { 01702 libmesh_assert(_sp.get()); 01703 libmesh_assert(!_n_nz || _n_nz == &_sp->n_nz); 01704 libmesh_assert(!_n_oz || _n_oz == &_sp->n_oz); 01705 _sp.reset(); 01706 } 01707 else 01708 { 01709 libmesh_assert(!_sp.get()); 01710 delete _n_nz; 01711 delete _n_oz; 01712 } 01713 _n_nz = NULL; 01714 _n_oz = NULL; 01715 } 01716 01717 01718 01719 void DofMap::extract_local_vector (const NumericVector<Number>& Ug, 01720 const std::vector<dof_id_type>& dof_indices_in, 01721 DenseVectorBase<Number>& Ue) const 01722 { 01723 #ifdef LIBMESH_ENABLE_AMR 01724 01725 // Trivial mapping 01726 libmesh_assert_equal_to (dof_indices_in.size(), Ue.size()); 01727 bool has_constrained_dofs = false; 01728 01729 for (unsigned int il=0; 01730 il != cast_int<unsigned int>(dof_indices_in.size()); il++) 01731 { 01732 const dof_id_type ig = dof_indices_in[il]; 01733 01734 if (this->is_constrained_dof (ig)) has_constrained_dofs = true; 01735 01736 libmesh_assert_less (ig, Ug.size()); 01737 01738 Ue.el(il) = Ug(ig); 01739 } 01740 01741 // If the element has any constrained DOFs then we need 01742 // to account for them in the mapping. This will handle 01743 // the case that the input vector is not constrained. 01744 if (has_constrained_dofs) 01745 { 01746 // Copy the input DOF indices. 01747 std::vector<dof_id_type> constrained_dof_indices(dof_indices_in); 01748 01749 DenseMatrix<Number> C; 01750 DenseVector<Number> H; 01751 01752 this->build_constraint_matrix_and_vector (C, H, constrained_dof_indices); 01753 01754 libmesh_assert_equal_to (dof_indices_in.size(), C.m()); 01755 libmesh_assert_equal_to (constrained_dof_indices.size(), C.n()); 01756 01757 // zero-out Ue 01758 Ue.zero(); 01759 01760 // compute Ue = C Ug, with proper mapping. 01761 const unsigned int n_original_dofs = 01762 cast_int<unsigned int>(dof_indices_in.size()); 01763 for (unsigned int i=0; i != n_original_dofs; i++) 01764 { 01765 Ue.el(i) = H(i); 01766 01767 const unsigned int n_constrained = 01768 cast_int<unsigned int>(constrained_dof_indices.size()); 01769 for (unsigned int j=0; j<n_constrained; j++) 01770 { 01771 const dof_id_type jg = constrained_dof_indices[j]; 01772 01773 // If Ug is a serial or ghosted vector, then this assert is 01774 // overzealous. If Ug is a parallel vector, then this assert 01775 // is redundant. 01776 // libmesh_assert ((jg >= Ug.first_local_index()) && 01777 // (jg < Ug.last_local_index())); 01778 01779 Ue.el(i) += C(i,j)*Ug(jg); 01780 } 01781 } 01782 } 01783 01784 #else 01785 01786 // Trivial mapping 01787 01788 const unsigned int n_original_dofs = 01789 cast_int<unsigned int>(dof_indices_in.size()); 01790 01791 libmesh_assert_equal_to (n_original_dofs, Ue.size()); 01792 01793 for (unsigned int il=0; il<n_original_dofs; il++) 01794 { 01795 const dof_id_type ig = dof_indices_in[il]; 01796 01797 libmesh_assert ((ig >= Ug.first_local_index()) && (ig < Ug.last_local_index())); 01798 01799 Ue.el(il) = Ug(ig); 01800 } 01801 01802 #endif 01803 } 01804 01805 void DofMap::dof_indices (const Elem* const elem, 01806 std::vector<dof_id_type>& di) const 01807 { 01808 // We now allow elem==NULL to request just SCALAR dofs 01809 // libmesh_assert(elem); 01810 01811 // If we are asking for current indices on an element, it ought to 01812 // be an active element (or a Side proxy, which also thinks it's 01813 // active) 01814 libmesh_assert(!elem || elem->active()); 01815 01816 START_LOG("dof_indices()", "DofMap"); 01817 01818 // Clear the DOF indices vector 01819 di.clear(); 01820 01821 const unsigned int n_vars = this->n_variables(); 01822 01823 #ifdef DEBUG 01824 // Check that sizes match in DEBUG mode 01825 std::size_t tot_size = 0; 01826 #endif 01827 01828 if (elem && elem->type() == TRI3SUBDIVISION) 01829 { 01830 // Subdivision surface FE require the 1-ring around elem 01831 const Tri3Subdivision* sd_elem = static_cast<const Tri3Subdivision*>(elem); 01832 01833 // Ghost subdivision elements have no real dofs 01834 if (!sd_elem->is_ghost()) 01835 { 01836 // Determine the nodes contributing to element elem 01837 std::vector<Node*> elem_nodes; 01838 MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes); 01839 01840 // Get the dof numbers 01841 for (unsigned int v=0; v<n_vars; v++) 01842 { 01843 if(this->variable(v).type().family == SCALAR && 01844 this->variable(v).active_on_subdomain(elem->subdomain_id())) 01845 { 01846 #ifdef DEBUG 01847 tot_size += this->variable(v).type().order; 01848 #endif 01849 std::vector<dof_id_type> di_new; 01850 this->SCALAR_dof_indices(di_new,v); 01851 di.insert( di.end(), di_new.begin(), di_new.end()); 01852 } 01853 else 01854 _dof_indices(elem, di, v, &elem_nodes[0], elem_nodes.size() 01855 #ifdef DEBUG 01856 , tot_size 01857 #endif 01858 ); 01859 } 01860 } 01861 01862 STOP_LOG("dof_indices()", "DofMap"); 01863 return; 01864 } 01865 01866 // Get the dof numbers 01867 for (unsigned int v=0; v<n_vars; v++) 01868 { 01869 const Variable & var = this->variable(v); 01870 if(var.type().family == SCALAR && 01871 (!elem || 01872 var.active_on_subdomain(elem->subdomain_id()))) 01873 { 01874 #ifdef DEBUG 01875 tot_size += var.type().order; 01876 #endif 01877 std::vector<dof_id_type> di_new; 01878 this->SCALAR_dof_indices(di_new,v); 01879 di.insert( di.end(), di_new.begin(), di_new.end()); 01880 } 01881 else if (elem) 01882 _dof_indices(elem, di, v, elem->get_nodes(), elem->n_nodes() 01883 #ifdef DEBUG 01884 , tot_size 01885 #endif 01886 ); 01887 } 01888 01889 #ifdef DEBUG 01890 libmesh_assert_equal_to (tot_size, di.size()); 01891 #endif 01892 01893 STOP_LOG("dof_indices()", "DofMap"); 01894 } 01895 01896 01897 void DofMap::dof_indices (const Elem* const elem, 01898 std::vector<dof_id_type>& di, 01899 const unsigned int vn) const 01900 { 01901 // We now allow elem==NULL to request just SCALAR dofs 01902 // libmesh_assert(elem); 01903 01904 START_LOG("dof_indices()", "DofMap"); 01905 01906 // Clear the DOF indices vector 01907 di.clear(); 01908 01909 #ifdef DEBUG 01910 // Check that sizes match in DEBUG mode 01911 std::size_t tot_size = 0; 01912 #endif 01913 01914 if (elem && elem->type() == TRI3SUBDIVISION) 01915 { 01916 // Subdivision surface FE require the 1-ring around elem 01917 const Tri3Subdivision* sd_elem = static_cast<const Tri3Subdivision*>(elem); 01918 01919 // Ghost subdivision elements have no real dofs 01920 if (!sd_elem->is_ghost()) 01921 { 01922 // Determine the nodes contributing to element elem 01923 std::vector<Node*> elem_nodes; 01924 MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes); 01925 01926 _dof_indices(elem, di, vn, &elem_nodes[0], elem_nodes.size() 01927 #ifdef DEBUG 01928 , tot_size 01929 #endif 01930 ); 01931 } 01932 01933 STOP_LOG("dof_indices()", "DofMap"); 01934 return; 01935 } 01936 01937 const Variable & var = this->variable(vn); 01938 01939 // Get the dof numbers 01940 if(var.type().family == SCALAR && 01941 (!elem || 01942 var.active_on_subdomain(elem->subdomain_id()))) 01943 { 01944 #ifdef DEBUG 01945 tot_size += var.type().order; 01946 #endif 01947 std::vector<dof_id_type> di_new; 01948 this->SCALAR_dof_indices(di_new,vn); 01949 di.insert( di.end(), di_new.begin(), di_new.end()); 01950 } 01951 else if (elem) 01952 _dof_indices(elem, di, vn, elem->get_nodes(), elem->n_nodes() 01953 #ifdef DEBUG 01954 , tot_size 01955 #endif 01956 ); 01957 01958 #ifdef DEBUG 01959 libmesh_assert_equal_to (tot_size, di.size()); 01960 #endif 01961 01962 STOP_LOG("dof_indices()", "DofMap"); 01963 } 01964 01965 01966 void DofMap::_dof_indices (const Elem* const elem, 01967 std::vector<dof_id_type>& di, 01968 const unsigned int v, 01969 const Node * const * nodes, 01970 unsigned int n_nodes 01971 #ifdef DEBUG 01972 , 01973 std::size_t & tot_size 01974 #endif 01975 ) const 01976 { 01977 // This internal function is only useful on valid elements 01978 libmesh_assert(elem); 01979 01980 const Variable & var = this->variable(v); 01981 01982 if (var.active_on_subdomain(elem->subdomain_id())) 01983 { 01984 const ElemType type = elem->type(); 01985 const unsigned int sys_num = this->sys_number(); 01986 const unsigned int dim = elem->dim(); 01987 01988 // Increase the polynomial order on p refined elements 01989 FEType fe_type = var.type(); 01990 fe_type.order = static_cast<Order>(fe_type.order + 01991 elem->p_level()); 01992 01993 const bool extra_hanging_dofs = 01994 FEInterface::extra_hanging_dofs(fe_type); 01995 01996 #ifdef DEBUG 01997 // The number of dofs per element is non-static for subdivision FE 01998 if (fe_type.family == SUBDIVISION) 01999 tot_size += n_nodes; 02000 else 02001 tot_size += FEInterface::n_dofs(dim,fe_type,type); 02002 #endif 02003 02004 // Get the node-based DOF numbers 02005 for (unsigned int n=0; n<n_nodes; n++) 02006 { 02007 const Node* node = nodes[n]; 02008 02009 // There is a potential problem with h refinement. Imagine a 02010 // quad9 that has a linear FE on it. Then, on the hanging side, 02011 // it can falsely identify a DOF at the mid-edge node. This is why 02012 // we call FEInterface instead of node->n_comp() directly. 02013 const unsigned int nc = FEInterface::n_dofs_at_node (dim, 02014 fe_type, 02015 type, 02016 n); 02017 02018 // If this is a non-vertex on a hanging node with extra 02019 // degrees of freedom, we use the non-vertex dofs (which 02020 // come in reverse order starting from the end, to 02021 // simplify p refinement) 02022 if (extra_hanging_dofs && !elem->is_vertex(n)) 02023 { 02024 const int dof_offset = node->n_comp(sys_num,v) - nc; 02025 02026 // We should never have fewer dofs than necessary on a 02027 // node unless we're getting indices on a parent element, 02028 // and we should never need the indices on such a node 02029 if (dof_offset < 0) 02030 { 02031 libmesh_assert(!elem->active()); 02032 di.resize(di.size() + nc, DofObject::invalid_id); 02033 } 02034 else 02035 for (int i=node->n_comp(sys_num,v)-1; i>=dof_offset; i--) 02036 { 02037 libmesh_assert_not_equal_to (node->dof_number(sys_num,v,i), 02038 DofObject::invalid_id); 02039 di.push_back(node->dof_number(sys_num,v,i)); 02040 } 02041 } 02042 // If this is a vertex or an element without extra hanging 02043 // dofs, our dofs come in forward order coming from the 02044 // beginning 02045 else 02046 for (unsigned int i=0; i<nc; i++) 02047 { 02048 libmesh_assert_not_equal_to (node->dof_number(sys_num,v,i), 02049 DofObject::invalid_id); 02050 di.push_back(node->dof_number(sys_num,v,i)); 02051 } 02052 } 02053 02054 // If there are any element-based DOF numbers, get them 02055 const unsigned int nc = FEInterface::n_dofs_per_elem(dim, 02056 fe_type, 02057 type); 02058 // We should never have fewer dofs than necessary on an 02059 // element unless we're getting indices on a parent element, 02060 // and we should never need those indices 02061 if (nc != 0) 02062 { 02063 if (elem->n_systems() > sys_num && 02064 nc <= elem->n_comp(sys_num,v)) 02065 { 02066 for (unsigned int i=0; i<nc; i++) 02067 { 02068 libmesh_assert_not_equal_to (elem->dof_number(sys_num,v,i), 02069 DofObject::invalid_id); 02070 02071 di.push_back(elem->dof_number(sys_num,v,i)); 02072 } 02073 } 02074 else 02075 { 02076 libmesh_assert(!elem->active() || fe_type.family == LAGRANGE || fe_type.family == SUBDIVISION); 02077 di.resize(di.size() + nc, DofObject::invalid_id); 02078 } 02079 } 02080 } 02081 } 02082 02083 02084 02085 void DofMap::SCALAR_dof_indices (std::vector<dof_id_type>& di, 02086 const unsigned int vn, 02087 #ifdef LIBMESH_ENABLE_AMR 02088 const bool old_dofs 02089 #else 02090 const bool 02091 #endif 02092 ) const 02093 { 02094 START_LOG("SCALAR_dof_indices()", "DofMap"); 02095 02096 libmesh_assert(this->variable(vn).type().family == SCALAR); 02097 02098 #ifdef LIBMESH_ENABLE_AMR 02099 // If we're asking for old dofs then we'd better have some 02100 if (old_dofs) 02101 libmesh_assert_greater_equal(n_old_dofs(), n_SCALAR_dofs()); 02102 02103 dof_id_type my_idx = old_dofs ? 02104 this->_first_old_scalar_df[vn] : this->_first_scalar_df[vn]; 02105 #else 02106 dof_id_type my_idx = this->_first_scalar_df[vn]; 02107 #endif 02108 02109 libmesh_assert_not_equal_to(my_idx, DofObject::invalid_id); 02110 02111 // The number of SCALAR dofs comes from the variable order 02112 const int n_dofs_vn = this->variable(vn).type().order; 02113 02114 di.resize(n_dofs_vn); 02115 for(int i = 0; i != n_dofs_vn; ++i) 02116 di[i] = my_idx++; 02117 02118 STOP_LOG("SCALAR_dof_indices()", "DofMap"); 02119 } 02120 02121 02122 02123 bool DofMap::all_semilocal_indices (const std::vector<dof_id_type>& dof_indices_in) const 02124 { 02125 // We're all semilocal unless we find a counterexample 02126 for (std::size_t i=0; i != dof_indices_in.size(); ++i) 02127 { 02128 const dof_id_type di = dof_indices_in[i]; 02129 // If it's not in the local indices 02130 if (di < this->first_dof() || 02131 di >= this->end_dof()) 02132 { 02133 // and if it's not in the ghost indices, then we're not 02134 // semilocal 02135 if (!std::binary_search(_send_list.begin(), _send_list.end(), di)) 02136 return false; 02137 } 02138 } 02139 02140 return true; 02141 } 02142 02143 02144 #ifdef LIBMESH_ENABLE_AMR 02145 02146 void DofMap::old_dof_indices (const Elem* const elem, 02147 std::vector<dof_id_type>& di, 02148 const unsigned int vn) const 02149 { 02150 START_LOG("old_dof_indices()", "DofMap"); 02151 02152 libmesh_assert(elem); 02153 02154 const ElemType type = elem->type(); 02155 const unsigned int sys_num = this->sys_number(); 02156 const unsigned int n_vars = this->n_variables(); 02157 const unsigned int dim = elem->dim(); 02158 02159 // If we have dof indices stored on the elem, and there's no chance 02160 // that we only have those indices because we were just p refined, 02161 // then we should have old dof indices too. 02162 libmesh_assert(!elem->has_dofs(sys_num) || 02163 elem->p_refinement_flag() == Elem::JUST_REFINED || 02164 elem->old_dof_object); 02165 02166 // Clear the DOF indices vector. 02167 di.clear(); 02168 02169 // Create a vector to indicate which 02170 // SCALAR variables have been requested 02171 std::vector<unsigned int> SCALAR_var_numbers; 02172 SCALAR_var_numbers.clear(); 02173 02174 // Determine the nodes contributing to element elem 02175 std::vector<Node*> elem_nodes; 02176 if (elem->type() == TRI3SUBDIVISION) 02177 { 02178 // Subdivision surface FE require the 1-ring around elem 02179 const Tri3Subdivision* sd_elem = static_cast<const Tri3Subdivision*>(elem); 02180 MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes); 02181 } 02182 else 02183 { 02184 // All other FE use only the nodes of elem itself 02185 elem_nodes.resize(elem->n_nodes(), NULL); 02186 for (unsigned int i=0; i<elem->n_nodes(); i++) 02187 elem_nodes[i] = elem->get_node(i); 02188 } 02189 02190 // Get the dof numbers 02191 for (unsigned int v=0; v<n_vars; v++) 02192 if ((v == vn) || (vn == libMesh::invalid_uint)) 02193 { 02194 if(this->variable(v).type().family == SCALAR && 02195 (!elem || 02196 this->variable(v).active_on_subdomain(elem->subdomain_id()))) 02197 { 02198 // We asked for this variable, so add it to the vector. 02199 std::vector<dof_id_type> di_new; 02200 this->SCALAR_dof_indices(di_new,v,true); 02201 di.insert( di.end(), di_new.begin(), di_new.end()); 02202 } 02203 else 02204 if (this->variable(v).active_on_subdomain(elem->subdomain_id())) 02205 { // Do this for all the variables if one was not specified 02206 // or just for the specified variable 02207 02208 // Increase the polynomial order on p refined elements, 02209 // but make sure you get the right polynomial order for 02210 // the OLD degrees of freedom 02211 int p_adjustment = 0; 02212 if (elem->p_refinement_flag() == Elem::JUST_REFINED) 02213 { 02214 libmesh_assert_greater (elem->p_level(), 0); 02215 p_adjustment = -1; 02216 } 02217 else if (elem->p_refinement_flag() == Elem::JUST_COARSENED) 02218 { 02219 p_adjustment = 1; 02220 } 02221 FEType fe_type = this->variable_type(v); 02222 fe_type.order = static_cast<Order>(fe_type.order + 02223 elem->p_level() + 02224 p_adjustment); 02225 02226 const bool extra_hanging_dofs = 02227 FEInterface::extra_hanging_dofs(fe_type); 02228 02229 // Get the node-based DOF numbers 02230 for (unsigned int n=0; n<elem_nodes.size(); n++) 02231 { 02232 const Node* node = elem_nodes[n]; 02233 02234 // There is a potential problem with h refinement. Imagine a 02235 // quad9 that has a linear FE on it. Then, on the hanging side, 02236 // it can falsely identify a DOF at the mid-edge node. This is why 02237 // we call FEInterface instead of node->n_comp() directly. 02238 const unsigned int nc = FEInterface::n_dofs_at_node (dim, 02239 fe_type, 02240 type, 02241 n); 02242 libmesh_assert(node->old_dof_object); 02243 02244 // If this is a non-vertex on a hanging node with extra 02245 // degrees of freedom, we use the non-vertex dofs (which 02246 // come in reverse order starting from the end, to 02247 // simplify p refinement) 02248 if (extra_hanging_dofs && !elem->is_vertex(n)) 02249 { 02250 const int dof_offset = 02251 node->old_dof_object->n_comp(sys_num,v) - nc; 02252 02253 // We should never have fewer dofs than necessary on a 02254 // node unless we're getting indices on a parent element 02255 // or a just-coarsened element 02256 if (dof_offset < 0) 02257 { 02258 libmesh_assert(!elem->active() || elem->refinement_flag() == 02259 Elem::JUST_COARSENED); 02260 di.resize(di.size() + nc, DofObject::invalid_id); 02261 } 02262 else 02263 for (int i=node->old_dof_object->n_comp(sys_num,v)-1; 02264 i>=dof_offset; i--) 02265 { 02266 libmesh_assert_not_equal_to (node->old_dof_object->dof_number(sys_num,v,i), 02267 DofObject::invalid_id); 02268 di.push_back(node->old_dof_object->dof_number(sys_num,v,i)); 02269 } 02270 } 02271 // If this is a vertex or an element without extra hanging 02272 // dofs, our dofs come in forward order coming from the 02273 // beginning 02274 else 02275 for (unsigned int i=0; i<nc; i++) 02276 { 02277 libmesh_assert_not_equal_to (node->old_dof_object->dof_number(sys_num,v,i), 02278 DofObject::invalid_id); 02279 di.push_back(node->old_dof_object->dof_number(sys_num,v,i)); 02280 } 02281 } 02282 02283 // If there are any element-based DOF numbers, get them 02284 const unsigned int nc = FEInterface::n_dofs_per_elem(dim, 02285 fe_type, 02286 type); 02287 02288 // We should never have fewer dofs than necessary on an 02289 // element unless we're getting indices on a parent element 02290 // or a just-coarsened element 02291 if (nc != 0) 02292 { 02293 if (elem->old_dof_object->n_systems() > sys_num && 02294 nc <= elem->old_dof_object->n_comp(sys_num,v)) 02295 { 02296 libmesh_assert(elem->old_dof_object); 02297 02298 for (unsigned int i=0; i<nc; i++) 02299 { 02300 libmesh_assert_not_equal_to (elem->old_dof_object->dof_number(sys_num,v,i), 02301 DofObject::invalid_id); 02302 02303 di.push_back(elem->old_dof_object->dof_number(sys_num,v,i)); 02304 } 02305 } 02306 else 02307 { 02308 libmesh_assert(!elem->active() || fe_type.family == LAGRANGE || 02309 elem->refinement_flag() == Elem::JUST_COARSENED); 02310 di.resize(di.size() + nc, DofObject::invalid_id); 02311 } 02312 } 02313 } 02314 } // end loop over variables 02315 02316 STOP_LOG("old_dof_indices()", "DofMap"); 02317 } 02318 02319 02320 02321 /* 02322 void DofMap::augment_send_list_for_projection(const MeshBase& mesh) 02323 { 02324 // Loop over the active local elements in the mesh. 02325 // If the element was just refined and its parent lives 02326 // on a different processor then we need to augment the 02327 // _send_list with the parent's DOF indices so that we 02328 // can access the parent's data for computing solution 02329 // projections, etc... 02330 02331 // The DOF indices for the parent 02332 std::vector<dof_id_type> di; 02333 02334 // Flag telling us if we need to re-sort the send_list. 02335 // Note we won't need to re-sort unless a child with a 02336 // parent belonging to a different processor is found 02337 bool needs_sorting = false; 02338 02339 02340 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 02341 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 02342 02343 for ( ; elem_it != elem_end; ++elem_it) 02344 { 02345 const Elem* const elem = *elem_it; 02346 02347 // We only need to consider the children that 02348 // were just refined 02349 if (elem->refinement_flag() != Elem::JUST_REFINED) continue; 02350 02351 const Elem* const parent = elem->parent(); 02352 02353 // If the parent lives on another processor 02354 // than the child 02355 if (parent != NULL) 02356 if (parent->processor_id() != elem->processor_id()) 02357 { 02358 // Get the DOF indices for the parent 02359 this->dof_indices (parent, di); 02360 02361 // Insert the DOF indices into the send list 02362 _send_list.insert (_send_list.end(), 02363 di.begin(), di.end()); 02364 02365 // We will need to re-sort the send list 02366 needs_sorting = true; 02367 } 02368 } 02369 02370 // The send-list might need to be sorted again. 02371 if (needs_sorting) this->sort_send_list (); 02372 } 02373 */ 02374 02375 #endif // LIBMESH_ENABLE_AMR 02376 02377 02378 #ifdef LIBMESH_ENABLE_CONSTRAINTS 02379 02380 void DofMap::find_connected_dofs (std::vector<dof_id_type>& elem_dofs) const 02381 { 02382 typedef std::set<dof_id_type> RCSet; 02383 02384 // First insert the DOFS we already depend on into the set. 02385 RCSet dof_set (elem_dofs.begin(), elem_dofs.end()); 02386 02387 bool done = true; 02388 02389 // Next insert any dofs those might be constrained in terms 02390 // of. Note that in this case we may not be done: Those may 02391 // in turn depend on others. So, we need to repeat this process 02392 // in that case until the system depends only on unconstrained 02393 // degrees of freedom. 02394 for (unsigned int i=0; i<elem_dofs.size(); i++) 02395 if (this->is_constrained_dof(elem_dofs[i])) 02396 { 02397 // If the DOF is constrained 02398 DofConstraints::const_iterator 02399 pos = _dof_constraints.find(elem_dofs[i]); 02400 02401 libmesh_assert (pos != _dof_constraints.end()); 02402 02403 const DofConstraintRow& constraint_row = pos->second; 02404 02405 // adaptive p refinement currently gives us lots of empty constraint 02406 // rows - we should optimize those DoFs away in the future. [RHS] 02407 //libmesh_assert (!constraint_row.empty()); 02408 02409 DofConstraintRow::const_iterator it = constraint_row.begin(); 02410 DofConstraintRow::const_iterator it_end = constraint_row.end(); 02411 02412 02413 // Add the DOFs this dof is constrained in terms of. 02414 // note that these dofs might also be constrained, so 02415 // we will need to call this function recursively. 02416 for ( ; it != it_end; ++it) 02417 if (!dof_set.count (it->first)) 02418 { 02419 dof_set.insert (it->first); 02420 done = false; 02421 } 02422 } 02423 02424 02425 // If not done then we need to do more work 02426 // (obviously :-) )! 02427 if (!done) 02428 { 02429 // Fill the vector with the contents of the set 02430 elem_dofs.clear(); 02431 elem_dofs.insert (elem_dofs.end(), 02432 dof_set.begin(), dof_set.end()); 02433 02434 02435 // May need to do this recursively. It is possible 02436 // that we just replaced a constrained DOF with another 02437 // constrained DOF. 02438 this->find_connected_dofs (elem_dofs); 02439 02440 } // end if (!done) 02441 } 02442 02443 #endif // LIBMESH_ENABLE_CONSTRAINTS 02444 02445 02446 02447 #if defined(__GNUC__) && (__GNUC__ < 4) && !defined(__INTEL_COMPILER) 02448 02449 void SparsityPattern::_dummy_function(void) 02450 { 02451 } 02452 02453 #endif 02454 02455 02456 02457 void SparsityPattern::Build::operator()(const ConstElemRange &range) 02458 { 02459 // Compute the sparsity structure of the global matrix. This can be 02460 // fed into a PetscMatrix to allocate exacly the number of nonzeros 02461 // necessary to store the matrix. This algorithm should be linear 02462 // in the (# of elements)*(# nodes per element) 02463 const processor_id_type proc_id = mesh.processor_id(); 02464 const dof_id_type n_dofs_on_proc = dof_map.n_dofs_on_processor(proc_id); 02465 const dof_id_type first_dof_on_proc = dof_map.first_dof(proc_id); 02466 const dof_id_type end_dof_on_proc = dof_map.end_dof(proc_id); 02467 02468 sparsity_pattern.resize(n_dofs_on_proc); 02469 02470 // If the user did not explicitly specify the DOF coupling 02471 // then all the DOFS are coupled to each other. Furthermore, 02472 // we can take a shortcut and do this more quickly here. So 02473 // we use an if-test. 02474 if ((dof_coupling == NULL) || (dof_coupling->empty())) 02475 { 02476 std::vector<dof_id_type> 02477 element_dofs, 02478 neighbor_dofs, 02479 dofs_to_add; 02480 02481 std::vector<const Elem*> active_neighbors; 02482 02483 for (ConstElemRange::const_iterator elem_it = range.begin() ; elem_it != range.end(); ++elem_it) 02484 { 02485 const Elem* const elem = *elem_it; 02486 02487 // Get the global indices of the DOFs with support on this element 02488 dof_map.dof_indices (elem, element_dofs); 02489 #ifdef LIBMESH_ENABLE_CONSTRAINTS 02490 dof_map.find_connected_dofs (element_dofs); 02491 #endif 02492 02493 // We can be more efficient if we sort the element DOFs 02494 // into increasing order 02495 std::sort(element_dofs.begin(), element_dofs.end()); 02496 02497 const unsigned int n_dofs_on_element = 02498 cast_int<unsigned int>(element_dofs.size()); 02499 02500 for (unsigned int i=0; i<n_dofs_on_element; i++) 02501 { 02502 const dof_id_type ig = element_dofs[i]; 02503 02504 SparsityPattern::Row *row; 02505 02506 // We save non-local row components for now so we can 02507 // communicate them to other processors later. 02508 02509 if ((ig >= first_dof_on_proc) && 02510 (ig < end_dof_on_proc)) 02511 { 02512 // This is what I mean 02513 // libmesh_assert_greater_equal ((ig - first_dof_on_proc), 0); 02514 // but do the test like this because ig and 02515 // first_dof_on_proc are unsigned ints 02516 libmesh_assert_greater_equal (ig, first_dof_on_proc); 02517 libmesh_assert_less ((ig - first_dof_on_proc), sparsity_pattern.size()); 02518 02519 row = &sparsity_pattern[ig - first_dof_on_proc]; 02520 } 02521 else 02522 { 02523 row = &nonlocal_pattern[ig]; 02524 } 02525 02526 // If the row is empty we will add *all* the element DOFs, 02527 // so just do that. 02528 if (row->empty()) 02529 { 02530 row->insert(row->end(), 02531 element_dofs.begin(), 02532 element_dofs.end()); 02533 } 02534 else 02535 { 02536 // Build a list of the DOF indices not found in the 02537 // sparsity pattern 02538 dofs_to_add.clear(); 02539 02540 // Cache iterators. Low will move forward, subsequent 02541 // searches will be on smaller ranges 02542 SparsityPattern::Row::iterator 02543 low = std::lower_bound (row->begin(), row->end(), element_dofs.front()), 02544 high = std::upper_bound (low, row->end(), element_dofs.back()); 02545 02546 for (unsigned int j=0; j<n_dofs_on_element; j++) 02547 { 02548 const dof_id_type jg = element_dofs[j]; 02549 02550 // See if jg is in the sorted range 02551 std::pair<SparsityPattern::Row::iterator, 02552 SparsityPattern::Row::iterator> 02553 pos = std::equal_range (low, high, jg); 02554 02555 // Must add jg if it wasn't found 02556 if (pos.first == pos.second) 02557 dofs_to_add.push_back(jg); 02558 02559 // pos.first is now a valid lower bound for any 02560 // remaining element DOFs. (That's why we sorted them.) 02561 // Use it for the next search 02562 low = pos.first; 02563 } 02564 02565 // Add to the sparsity pattern 02566 if (!dofs_to_add.empty()) 02567 { 02568 const std::size_t old_size = row->size(); 02569 02570 row->insert (row->end(), 02571 dofs_to_add.begin(), 02572 dofs_to_add.end()); 02573 02574 SparsityPattern::sort_row 02575 (row->begin(), row->begin()+old_size, row->end()); 02576 } 02577 } 02578 02579 // Now (possibly) add dofs from neighboring elements 02580 // TODO:[BSK] optimize this like above! 02581 if (implicit_neighbor_dofs) 02582 for (unsigned int s=0; s<elem->n_sides(); s++) 02583 if (elem->neighbor(s) != NULL) 02584 { 02585 const Elem* const neighbor_0 = elem->neighbor(s); 02586 #ifdef LIBMESH_ENABLE_AMR 02587 neighbor_0->active_family_tree_by_neighbor(active_neighbors,elem); 02588 #else 02589 active_neighbors.clear(); 02590 active_neighbors.push_back(neighbor_0); 02591 #endif 02592 02593 for (std::size_t a=0; a != active_neighbors.size(); ++a) 02594 { 02595 const Elem *neighbor = active_neighbors[a]; 02596 02597 dof_map.dof_indices (neighbor, neighbor_dofs); 02598 #ifdef LIBMESH_ENABLE_CONSTRAINTS 02599 dof_map.find_connected_dofs (neighbor_dofs); 02600 #endif 02601 const std::size_t n_dofs_on_neighbor = neighbor_dofs.size(); 02602 02603 for (std::size_t j=0; j<n_dofs_on_neighbor; j++) 02604 { 02605 const dof_id_type jg = neighbor_dofs[j]; 02606 02607 // See if jg is in the sorted range 02608 std::pair<SparsityPattern::Row::iterator, 02609 SparsityPattern::Row::iterator> 02610 pos = std::equal_range (row->begin(), row->end(), jg); 02611 02612 // Insert jg if it wasn't found 02613 if (pos.first == pos.second) 02614 row->insert (pos.first, jg); 02615 } 02616 } 02617 } 02618 } 02619 } 02620 } 02621 02622 // This is what we do in the case that the user has specified 02623 // explicit DOF coupling. 02624 else 02625 { 02626 libmesh_assert(dof_coupling); 02627 libmesh_assert_equal_to (dof_coupling->size(), 02628 dof_map.n_variables()); 02629 02630 const unsigned int n_var = dof_map.n_variables(); 02631 02632 std::vector<dof_id_type> 02633 element_dofs_i, 02634 element_dofs_j, 02635 neighbor_dofs, 02636 dofs_to_add; 02637 02638 02639 std::vector<const Elem*> active_neighbors; 02640 for (ConstElemRange::const_iterator elem_it = range.begin() ; elem_it != range.end(); ++elem_it) 02641 for (unsigned int vi=0; vi<n_var; vi++) 02642 { 02643 const Elem* const elem = *elem_it; 02644 02645 // Find element dofs for variable vi 02646 dof_map.dof_indices (elem, element_dofs_i, vi); 02647 #ifdef LIBMESH_ENABLE_CONSTRAINTS 02648 dof_map.find_connected_dofs (element_dofs_i); 02649 #endif 02650 02651 // We can be more efficient if we sort the element DOFs 02652 // into increasing order 02653 std::sort(element_dofs_i.begin(), element_dofs_i.end()); 02654 const unsigned int n_dofs_on_element_i = 02655 cast_int<unsigned int>(element_dofs_i.size()); 02656 02657 for (unsigned int vj=0; vj<n_var; vj++) 02658 if ((*dof_coupling)(vi,vj)) // If vi couples to vj 02659 { 02660 // Find element dofs for variable vj, note that 02661 // if vi==vj we already have the dofs. 02662 if (vi != vj) 02663 { 02664 dof_map.dof_indices (elem, element_dofs_j, vj); 02665 #ifdef LIBMESH_ENABLE_CONSTRAINTS 02666 dof_map.find_connected_dofs (element_dofs_j); 02667 #endif 02668 02669 // We can be more efficient if we sort the element DOFs 02670 // into increasing order 02671 std::sort (element_dofs_j.begin(), element_dofs_j.end()); 02672 } 02673 else 02674 element_dofs_j = element_dofs_i; 02675 02676 const unsigned int n_dofs_on_element_j = 02677 cast_int<unsigned int>(element_dofs_j.size()); 02678 02679 // there might be 0 dofs for the other variable on the same element (when subdomain variables do not overlap) and that's when we do not do anything 02680 if (n_dofs_on_element_j > 0) 02681 { 02682 for (unsigned int i=0; i<n_dofs_on_element_i; i++) 02683 { 02684 const dof_id_type ig = element_dofs_i[i]; 02685 02686 SparsityPattern::Row *row; 02687 02688 // We save non-local row components for now so we can 02689 // communicate them to other processors later. 02690 02691 if ((ig >= first_dof_on_proc) && 02692 (ig < end_dof_on_proc)) 02693 { 02694 // This is what I mean 02695 // libmesh_assert_greater_equal ((ig - first_dof_on_proc), 0); 02696 // but do the test like this because ig and 02697 // first_dof_on_proc are unsigned ints 02698 libmesh_assert_greater_equal (ig, first_dof_on_proc); 02699 libmesh_assert_less (ig, (sparsity_pattern.size() + 02700 first_dof_on_proc)); 02701 02702 row = &sparsity_pattern[ig - first_dof_on_proc]; 02703 } 02704 else 02705 { 02706 row = &nonlocal_pattern[ig]; 02707 } 02708 02709 // If the row is empty we will add *all* the element j DOFs, 02710 // so just do that. 02711 if (row->empty()) 02712 { 02713 row->insert(row->end(), 02714 element_dofs_j.begin(), 02715 element_dofs_j.end()); 02716 } 02717 else 02718 { 02719 // Build a list of the DOF indices not found in the 02720 // sparsity pattern 02721 dofs_to_add.clear(); 02722 02723 // Cache iterators. Low will move forward, subsequent 02724 // searches will be on smaller ranges 02725 SparsityPattern::Row::iterator 02726 low = std::lower_bound 02727 (row->begin(), row->end(), element_dofs_j.front()), 02728 high = std::upper_bound 02729 (low, row->end(), element_dofs_j.back()); 02730 02731 for (unsigned int j=0; j<n_dofs_on_element_j; j++) 02732 { 02733 const dof_id_type jg = element_dofs_j[j]; 02734 02735 // See if jg is in the sorted range 02736 std::pair<SparsityPattern::Row::iterator, 02737 SparsityPattern::Row::iterator> 02738 pos = std::equal_range (low, high, jg); 02739 02740 // Must add jg if it wasn't found 02741 if (pos.first == pos.second) 02742 dofs_to_add.push_back(jg); 02743 02744 // pos.first is now a valid lower bound for any 02745 // remaining element j DOFs. (That's why we sorted them.) 02746 // Use it for the next search 02747 low = pos.first; 02748 } 02749 02750 // Add to the sparsity pattern 02751 if (!dofs_to_add.empty()) 02752 { 02753 const std::size_t old_size = row->size(); 02754 02755 row->insert (row->end(), 02756 dofs_to_add.begin(), 02757 dofs_to_add.end()); 02758 02759 SparsityPattern::sort_row 02760 (row->begin(), row->begin()+old_size, 02761 row->end()); 02762 } 02763 } 02764 // Now (possibly) add dofs from neighboring elements 02765 // TODO:[BSK] optimize this like above! 02766 if (implicit_neighbor_dofs) 02767 for (unsigned int s=0; s<elem->n_sides(); s++) 02768 if (elem->neighbor(s) != NULL) 02769 { 02770 const Elem* const neighbor_0 = elem->neighbor(s); 02771 #ifdef LIBMESH_ENABLE_AMR 02772 neighbor_0->active_family_tree_by_neighbor(active_neighbors,elem); 02773 #else 02774 active_neighbors.clear(); 02775 active_neighbors.push_back(neighbor_0); 02776 #endif 02777 02778 for (std::size_t a=0; a != active_neighbors.size(); ++a) 02779 { 02780 const Elem *neighbor = active_neighbors[a]; 02781 02782 dof_map.dof_indices (neighbor, neighbor_dofs); 02783 #ifdef LIBMESH_ENABLE_CONSTRAINTS 02784 dof_map.find_connected_dofs (neighbor_dofs); 02785 #endif 02786 const unsigned int n_dofs_on_neighbor = 02787 cast_int<unsigned int>(neighbor_dofs.size()); 02788 02789 for (unsigned int j=0; j<n_dofs_on_neighbor; j++) 02790 { 02791 const dof_id_type jg = neighbor_dofs[j]; 02792 02793 // See if jg is in the sorted range 02794 std::pair<SparsityPattern::Row::iterator, 02795 SparsityPattern::Row::iterator> 02796 pos = std::equal_range (row->begin(), row->end(), jg); 02797 02798 // Insert jg if it wasn't found 02799 if (pos.first == pos.second) 02800 row->insert (pos.first, jg); 02801 } 02802 } 02803 } 02804 } 02805 } 02806 } 02807 } 02808 } 02809 02810 // Now a new chunk of sparsity structure is built for all of the 02811 // DOFs connected to our rows of the matrix. 02812 02813 // If we're building a full sparsity pattern, then we've got 02814 // complete rows to work with, so we can just count them from 02815 // scratch. 02816 if (need_full_sparsity_pattern) 02817 { 02818 n_nz.clear(); 02819 n_oz.clear(); 02820 } 02821 02822 n_nz.resize (n_dofs_on_proc, 0); 02823 n_oz.resize (n_dofs_on_proc, 0); 02824 02825 for (dof_id_type i=0; i<n_dofs_on_proc; i++) 02826 { 02827 // Get the row of the sparsity pattern 02828 SparsityPattern::Row &row = sparsity_pattern[i]; 02829 02830 for (dof_id_type j=0; j<row.size(); j++) 02831 if ((row[j] < first_dof_on_proc) || (row[j] >= end_dof_on_proc)) 02832 n_oz[i]++; 02833 else 02834 n_nz[i]++; 02835 02836 // If we're not building a full sparsity pattern, then we want 02837 // to avoid overcounting these entries as much as possible. 02838 if (!need_full_sparsity_pattern) 02839 row.clear(); 02840 } 02841 } 02842 02843 02844 02845 void SparsityPattern::Build::join (const SparsityPattern::Build &other) 02846 { 02847 const processor_id_type proc_id = mesh.processor_id(); 02848 const dof_id_type n_global_dofs = dof_map.n_dofs(); 02849 const dof_id_type n_dofs_on_proc = dof_map.n_dofs_on_processor(proc_id); 02850 const dof_id_type first_dof_on_proc = dof_map.first_dof(proc_id); 02851 const dof_id_type end_dof_on_proc = dof_map.end_dof(proc_id); 02852 02853 libmesh_assert_equal_to (sparsity_pattern.size(), other.sparsity_pattern.size()); 02854 libmesh_assert_equal_to (n_nz.size(), sparsity_pattern.size()); 02855 libmesh_assert_equal_to (n_oz.size(), sparsity_pattern.size()); 02856 02857 for (dof_id_type r=0; r<n_dofs_on_proc; r++) 02858 { 02859 // increment the number of on and off-processor nonzeros in this row 02860 // (note this will be an upper bound unless we need the full sparsity pattern) 02861 if (need_full_sparsity_pattern) 02862 { 02863 SparsityPattern::Row &my_row = sparsity_pattern[r]; 02864 const SparsityPattern::Row &their_row = other.sparsity_pattern[r]; 02865 02866 // simple copy if I have no dofs 02867 if (my_row.empty()) 02868 my_row = their_row; 02869 02870 // otherwise add their DOFs to mine, resort, and re-unique the row 02871 else if (!their_row.empty()) // do nothing for the trivial case where 02872 { // their row is empty 02873 my_row.insert (my_row.end(), 02874 their_row.begin(), 02875 their_row.end()); 02876 02877 // We cannot use SparsityPattern::sort_row() here because it expects 02878 // the [begin,middle) [middle,end) to be non-overlapping. This is not 02879 // necessarily the case here, so use std::sort() 02880 std::sort (my_row.begin(), my_row.end()); 02881 02882 my_row.erase(std::unique (my_row.begin(), my_row.end()), my_row.end()); 02883 } 02884 02885 // fix the number of on and off-processor nonzeros in this row 02886 n_nz[r] = n_oz[r] = 0; 02887 02888 for (std::size_t j=0; j<my_row.size(); j++) 02889 if ((my_row[j] < first_dof_on_proc) || (my_row[j] >= end_dof_on_proc)) 02890 n_oz[r]++; 02891 else 02892 n_nz[r]++; 02893 } 02894 else 02895 { 02896 n_nz[r] += other.n_nz[r]; 02897 n_nz[r] = std::min(n_nz[r], n_dofs_on_proc); 02898 n_oz[r] += other.n_oz[r]; 02899 n_oz[r] =std::min(n_oz[r], static_cast<dof_id_type>(n_global_dofs-n_nz[r])); 02900 } 02901 } 02902 02903 // Move nonlocal row information to ourselves; the other thread 02904 // won't need it in the map after that. 02905 NonlocalGraph::const_iterator it = other.nonlocal_pattern.begin(); 02906 for (; it != other.nonlocal_pattern.end(); ++it) 02907 { 02908 #ifndef NDEBUG 02909 const dof_id_type dof_id = it->first; 02910 02911 processor_id_type dbg_proc_id = 0; 02912 while (dof_id >= dof_map.end_dof(dbg_proc_id)) 02913 dbg_proc_id++; 02914 libmesh_assert (dbg_proc_id != this->processor_id()); 02915 #endif 02916 02917 const SparsityPattern::Row &their_row = it->second; 02918 02919 // We should have no empty values in a map 02920 libmesh_assert (!their_row.empty()); 02921 02922 NonlocalGraph::iterator my_it = nonlocal_pattern.find(it->first); 02923 if (my_it == nonlocal_pattern.end()) 02924 { 02925 // nonlocal_pattern[it->first].swap(their_row); 02926 nonlocal_pattern[it->first] = their_row; 02927 } 02928 else 02929 { 02930 SparsityPattern::Row &my_row = my_it->second; 02931 02932 my_row.insert (my_row.end(), 02933 their_row.begin(), 02934 their_row.end()); 02935 02936 // We cannot use SparsityPattern::sort_row() here because it expects 02937 // the [begin,middle) [middle,end) to be non-overlapping. This is not 02938 // necessarily the case here, so use std::sort() 02939 std::sort (my_row.begin(), my_row.end()); 02940 02941 my_row.erase(std::unique (my_row.begin(), my_row.end()), my_row.end()); 02942 } 02943 } 02944 } 02945 02946 02947 02948 void SparsityPattern::Build::parallel_sync () 02949 { 02950 parallel_object_only(); 02951 this->comm().verify(need_full_sparsity_pattern); 02952 02953 const dof_id_type n_global_dofs = dof_map.n_dofs(); 02954 const dof_id_type n_dofs_on_proc = dof_map.n_dofs_on_processor(this->processor_id()); 02955 const dof_id_type local_first_dof = dof_map.first_dof(); 02956 const dof_id_type local_end_dof = dof_map.end_dof(); 02957 02958 // Trade sparsity rows with other processors 02959 for (processor_id_type p=1; p != this->n_processors(); ++p) 02960 { 02961 // Push to processor procup while receiving from procdown 02962 processor_id_type procup = 02963 cast_int<processor_id_type>((this->processor_id() + p) % 02964 this->n_processors()); 02965 processor_id_type procdown = 02966 cast_int<processor_id_type>((this->n_processors() + this->processor_id() - p) % 02967 this->n_processors()); 02968 02969 // Pack the sparsity pattern rows to push to procup 02970 std::vector<dof_id_type> pushed_row_ids, 02971 pushed_row_ids_to_me; 02972 std::vector<std::vector<dof_id_type> > pushed_rows, 02973 pushed_rows_to_me; 02974 02975 // Move nonlocal row information to a structure to send it from; 02976 // we don't need it in the map after that. 02977 NonlocalGraph::iterator it = nonlocal_pattern.begin(); 02978 while (it != nonlocal_pattern.end()) 02979 { 02980 const dof_id_type dof_id = it->first; 02981 processor_id_type proc_id = 0; 02982 while (dof_id >= dof_map.end_dof(proc_id)) 02983 proc_id++; 02984 02985 libmesh_assert (proc_id != this->processor_id()); 02986 02987 if (proc_id == procup) 02988 { 02989 pushed_row_ids.push_back(dof_id); 02990 02991 // We can't just do the swap trick here, thanks to the 02992 // differing vector allocators? 02993 pushed_rows.push_back(std::vector<dof_id_type>()); 02994 pushed_rows.back().assign 02995 (it->second.begin(), it->second.end()); 02996 02997 nonlocal_pattern.erase(it++); 02998 } 02999 else 03000 ++it; 03001 } 03002 03003 this->comm().send_receive(procup, pushed_row_ids, 03004 procdown, pushed_row_ids_to_me); 03005 this->comm().send_receive(procup, pushed_rows, 03006 procdown, pushed_rows_to_me); 03007 pushed_row_ids.clear(); 03008 pushed_rows.clear(); 03009 03010 const std::size_t n_rows = pushed_row_ids_to_me.size(); 03011 for (std::size_t i=0; i != n_rows; ++i) 03012 { 03013 const dof_id_type r = pushed_row_ids_to_me[i]; 03014 const dof_id_type my_r = r - local_first_dof; 03015 03016 std::vector<dof_id_type> &their_row = pushed_rows_to_me[i]; 03017 03018 if (need_full_sparsity_pattern) 03019 { 03020 SparsityPattern::Row &my_row = 03021 sparsity_pattern[my_r]; 03022 03023 // They wouldn't have sent an empty row 03024 libmesh_assert(!their_row.empty()); 03025 03026 // We can end up with an empty row on a dof that touches our 03027 // inactive elements but not our active ones 03028 if (my_row.empty()) 03029 { 03030 my_row.assign (their_row.begin(), 03031 their_row.end()); 03032 } 03033 else 03034 { 03035 my_row.insert (my_row.end(), 03036 their_row.begin(), 03037 their_row.end()); 03038 03039 // We cannot use SparsityPattern::sort_row() here because it expects 03040 // the [begin,middle) [middle,end) to be non-overlapping. This is not 03041 // necessarily the case here, so use std::sort() 03042 std::sort (my_row.begin(), my_row.end()); 03043 03044 my_row.erase(std::unique (my_row.begin(), my_row.end()), my_row.end()); 03045 } 03046 03047 // fix the number of on and off-processor nonzeros in this row 03048 n_nz[my_r] = n_oz[my_r] = 0; 03049 03050 for (std::size_t j=0; j<my_row.size(); j++) 03051 if ((my_row[j] < local_first_dof) || (my_row[j] >= local_end_dof)) 03052 n_oz[my_r]++; 03053 else 03054 n_nz[my_r]++; 03055 } 03056 else 03057 { 03058 for (std::size_t j=0; j<their_row.size(); j++) 03059 if ((their_row[j] < local_first_dof) || (their_row[j] >= local_end_dof)) 03060 n_oz[my_r]++; 03061 else 03062 n_nz[my_r]++; 03063 03064 n_nz[my_r] = std::min(n_nz[my_r], n_dofs_on_proc); 03065 n_oz[my_r] = std::min(n_oz[my_r], 03066 static_cast<dof_id_type>(n_global_dofs-n_nz[my_r])); 03067 } 03068 } 03069 } 03070 03071 // We should have sent everything at this point. 03072 libmesh_assert (nonlocal_pattern.empty()); 03073 } 03074 03075 03076 03077 void DofMap::print_info(std::ostream& os) const 03078 { 03079 os << this->get_info(); 03080 } 03081 03082 03083 03084 std::string DofMap::get_info() const 03085 { 03086 std::ostringstream os; 03087 03088 // If we didn't calculate the exact sparsity pattern, the threaded 03089 // sparsity pattern assembly may have just given us an upper bound 03090 // on sparsity. 03091 const char* may_equal = " <= "; 03092 03093 // If we calculated the exact sparsity pattern, then we can report 03094 // exact bandwidth figures: 03095 std::vector<SparseMatrix<Number>* >::const_iterator 03096 pos = _matrices.begin(), 03097 end = _matrices.end(); 03098 03099 for (; pos != end; ++pos) 03100 if ((*pos)->need_full_sparsity_pattern()) 03101 may_equal = " = "; 03102 03103 dof_id_type max_n_nz = 0, max_n_oz = 0; 03104 long double avg_n_nz = 0, avg_n_oz = 0; 03105 03106 if (_n_nz) 03107 { 03108 for (std::size_t i = 0; i != _n_nz->size(); ++i) 03109 { 03110 max_n_nz = std::max(max_n_nz, (*_n_nz)[i]); 03111 avg_n_nz += (*_n_nz)[i]; 03112 } 03113 03114 std::size_t n_nz_size = _n_nz->size(); 03115 03116 this->comm().max(max_n_nz); 03117 this->comm().sum(avg_n_nz); 03118 this->comm().sum(n_nz_size); 03119 03120 avg_n_nz /= std::max(n_nz_size,std::size_t(1)); 03121 03122 libmesh_assert(_n_oz); 03123 03124 for (std::size_t i = 0; i != (*_n_oz).size(); ++i) 03125 { 03126 max_n_oz = std::max(max_n_oz, (*_n_oz)[i]); 03127 avg_n_oz += (*_n_oz)[i]; 03128 } 03129 03130 std::size_t n_oz_size = _n_oz->size(); 03131 03132 this->comm().max(max_n_oz); 03133 this->comm().sum(avg_n_oz); 03134 this->comm().sum(n_oz_size); 03135 03136 avg_n_oz /= std::max(n_oz_size,std::size_t(1)); 03137 } 03138 03139 os << " DofMap Sparsity\n Average On-Processor Bandwidth" 03140 << may_equal << avg_n_nz << '\n'; 03141 03142 os << " Average Off-Processor Bandwidth" 03143 << may_equal << avg_n_oz << '\n'; 03144 03145 os << " Maximum On-Processor Bandwidth" 03146 << may_equal << max_n_nz << '\n'; 03147 03148 os << " Maximum Off-Processor Bandwidth" 03149 << may_equal << max_n_oz << std::endl; 03150 03151 #ifdef LIBMESH_ENABLE_CONSTRAINTS 03152 03153 std::size_t n_constraints = 0, max_constraint_length = 0, 03154 n_rhss = 0; 03155 long double avg_constraint_length = 0.; 03156 03157 for (DofConstraints::const_iterator it=_dof_constraints.begin(); 03158 it != _dof_constraints.end(); ++it) 03159 { 03160 // Only count local constraints, then sum later 03161 const dof_id_type constrained_dof = it->first; 03162 if (constrained_dof < this->first_dof() || 03163 constrained_dof >= this->end_dof()) 03164 continue; 03165 03166 const DofConstraintRow& row = it->second; 03167 std::size_t rowsize = row.size(); 03168 03169 max_constraint_length = std::max(max_constraint_length, 03170 rowsize); 03171 avg_constraint_length += rowsize; 03172 n_constraints++; 03173 03174 if (_primal_constraint_values.count(constrained_dof)) 03175 n_rhss++; 03176 } 03177 03178 this->comm().sum(n_constraints); 03179 this->comm().sum(n_rhss); 03180 this->comm().sum(avg_constraint_length); 03181 this->comm().max(max_constraint_length); 03182 03183 os << " DofMap Constraints\n Number of DoF Constraints = " 03184 << n_constraints; 03185 if (n_rhss) 03186 os << '\n' 03187 << " Number of Heterogenous Constraints= " << n_rhss; 03188 if (n_constraints) 03189 { 03190 avg_constraint_length /= n_constraints; 03191 03192 os << '\n' 03193 << " Average DoF Constraint Length= " << avg_constraint_length; 03194 } 03195 03196 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 03197 std::size_t n_node_constraints = 0, max_node_constraint_length = 0, 03198 n_node_rhss = 0; 03199 long double avg_node_constraint_length = 0.; 03200 03201 for (NodeConstraints::const_iterator it=_node_constraints.begin(); 03202 it != _node_constraints.end(); ++it) 03203 { 03204 // Only count local constraints, then sum later 03205 const Node *node = it->first; 03206 if (node->processor_id() != this->processor_id()) 03207 continue; 03208 03209 const NodeConstraintRow& row = it->second.first; 03210 std::size_t rowsize = row.size(); 03211 03212 max_node_constraint_length = std::max(max_node_constraint_length, 03213 rowsize); 03214 avg_node_constraint_length += rowsize; 03215 n_node_constraints++; 03216 03217 if (it->second.second != Point(0)) 03218 n_node_rhss++; 03219 } 03220 03221 this->comm().sum(n_node_constraints); 03222 this->comm().sum(n_node_rhss); 03223 this->comm().sum(avg_node_constraint_length); 03224 this->comm().max(max_node_constraint_length); 03225 03226 os << "\n Number of Node Constraints = " << n_node_constraints; 03227 if (n_node_rhss) 03228 os << '\n' 03229 << " Number of Heterogenous Node Constraints= " << n_node_rhss; 03230 if (n_node_constraints) 03231 { 03232 avg_node_constraint_length /= n_node_constraints; 03233 os << "\n Maximum Node Constraint Length= " << max_node_constraint_length 03234 << '\n' 03235 << " Average Node Constraint Length= " << avg_node_constraint_length; 03236 } 03237 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 03238 03239 os << std::endl; 03240 03241 #endif // LIBMESH_ENABLE_CONSTRAINTS 03242 03243 return os.str(); 03244 } 03245 03246 } // namespace libMesh