$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 #ifndef LIBMESH_DOF_OBJECT_H 00021 #define LIBMESH_DOF_OBJECT_H 00022 00023 // Local includes 00024 #include "libmesh/id_types.h" 00025 #include "libmesh/libmesh_config.h" 00026 #include "libmesh/libmesh_common.h" 00027 #include "libmesh/libmesh.h" // libMesh::invalid_uint 00028 #include "libmesh/reference_counted_object.h" 00029 00030 // C++ includes 00031 #include <cstddef> 00032 #include <vector> 00033 00034 namespace libMesh 00035 { 00036 00037 // Forward declarations 00038 class DofObject; 00039 00040 00051 class DofObject : public ReferenceCountedObject<DofObject> 00052 { 00053 #ifdef LIBMESH_IS_UNIT_TESTING 00054 public: 00055 #else 00056 protected: 00057 #endif 00058 00063 DofObject (); 00064 00069 ~DofObject (); 00070 00074 DofObject (const DofObject&); 00075 00079 DofObject& operator= (const DofObject& dof_obj); 00080 00081 public: 00082 00083 #ifdef LIBMESH_ENABLE_AMR 00084 00089 DofObject* old_dof_object; 00090 00094 void clear_old_dof_object (); 00095 00099 void set_old_dof_object (); 00100 00101 #endif 00102 00107 void clear_dofs (); 00108 00112 void invalidate_dofs (const unsigned int sys_num = libMesh::invalid_uint); 00113 00117 void invalidate_id (); 00118 00122 void invalidate_processor_id (); 00123 00127 void invalidate (); 00128 00134 unsigned int n_dofs (const unsigned int s, 00135 const unsigned int var = 00136 libMesh::invalid_uint) const; 00137 00141 dof_id_type id () const; 00142 00146 dof_id_type & set_id (); 00147 00151 unique_id_type unique_id () const; 00152 00156 unique_id_type & set_unique_id (); 00157 00161 void set_id (const dof_id_type dofid) 00162 { this->set_id() = dofid; } 00163 00168 bool valid_id () const; 00169 00174 bool valid_unique_id () const; 00175 00183 processor_id_type processor_id () const; 00184 00189 processor_id_type& processor_id (); 00190 00194 void processor_id (const processor_id_type pid); 00195 00200 bool valid_processor_id () const; 00201 00206 unsigned int n_systems() const; 00207 00211 void set_n_systems (const unsigned int s); 00212 00216 void add_system (); 00217 00222 unsigned int n_var_groups(const unsigned int s) const; 00223 00228 unsigned int n_vars(const unsigned int s, 00229 const unsigned int vg) const; 00230 00235 unsigned int n_vars(const unsigned int s) const; 00236 00244 void set_n_vars_per_group(const unsigned int s, 00245 const std::vector<unsigned int> &nvpg); 00246 00256 unsigned int n_comp(const unsigned int s, 00257 const unsigned int var) const; 00258 00268 unsigned int n_comp_group(const unsigned int s, 00269 const unsigned int vg) const; 00270 00275 void set_n_comp(const unsigned int s, 00276 const unsigned int var, 00277 const unsigned int ncomp); 00278 00283 void set_n_comp_group(const unsigned int s, 00284 const unsigned int vg, 00285 const unsigned int ncomp); 00286 00296 dof_id_type dof_number(const unsigned int s, 00297 const unsigned int var, 00298 const unsigned int comp) const; 00299 00304 void set_dof_number(const unsigned int s, 00305 const unsigned int var, 00306 const unsigned int comp, 00307 const dof_id_type dn); 00308 00313 bool has_dofs(const unsigned int s=libMesh::invalid_uint) const; 00314 00320 void set_vg_dof_base(const unsigned int s, 00321 const unsigned int vg, 00322 const dof_id_type db); 00323 00329 dof_id_type vg_dof_base(const unsigned int s, 00330 const unsigned int vg) const; 00331 00335 static const dof_id_type invalid_id = static_cast<dof_id_type>(-1); 00336 00340 static const unique_id_type invalid_unique_id = static_cast<unique_id_type>(-1); 00341 00346 static const processor_id_type invalid_processor_id = static_cast<processor_id_type>(-1); 00347 00352 unsigned int packed_indexing_size() const; 00353 00358 static unsigned int unpackable_indexing_size 00359 (std::vector<largest_id_type>::const_iterator begin); 00360 00366 void unpack_indexing(std::vector<largest_id_type>::const_iterator begin); 00367 00373 void pack_indexing(std::back_insert_iterator<std::vector<largest_id_type> > target) const; 00374 00378 void debug_buffer () const; 00379 00380 private: 00381 00386 unsigned int var_to_vg (const unsigned int s, 00387 const unsigned int var) const; 00388 00393 unsigned int system_var_to_vg_var (const unsigned int s, 00394 const unsigned int vg, 00395 const unsigned int var) const; 00396 00400 #ifdef LIBMESH_ENABLE_UNIQUE_ID 00401 unique_id_type _unique_id; 00402 #endif 00403 00407 dof_id_type _id; 00408 00418 processor_id_type _processor_id; 00419 00470 typedef dof_id_type index_t; 00471 typedef std::vector<index_t> index_buffer_t; 00472 index_buffer_t _idx_buf; 00473 00483 static const index_t ncv_magic = 256; // = 2^8, in case we want to manually bitshift 00484 static const index_t ncv_magic_exp = 8; // Let's manually bitshift 00485 00489 unsigned int start_idx(const unsigned int s) const; 00490 00494 unsigned int end_idx(const unsigned int s) const; 00495 00496 // methods only available for unit testing 00497 #ifdef LIBMESH_IS_UNIT_TESTING 00498 public: 00499 void set_buffer (const std::vector<dof_id_type> &buf) 00500 { _idx_buf = buf; } 00501 #endif 00502 }; 00503 00504 00505 00506 //------------------------------------------------------ 00507 // Inline functions 00508 inline 00509 DofObject::DofObject () : 00510 #ifdef LIBMESH_ENABLE_AMR 00511 old_dof_object(NULL), 00512 #endif 00513 #ifdef LIBMESH_ENABLE_UNIQUE_ID 00514 _unique_id (invalid_unique_id), 00515 #endif 00516 _id (invalid_id), 00517 _processor_id (invalid_processor_id) 00518 { 00519 this->invalidate(); 00520 } 00521 00522 00523 00524 00525 00526 inline 00527 DofObject::~DofObject () 00528 { 00529 // Free all memory. 00530 #ifdef LIBMESH_ENABLE_AMR 00531 this->clear_old_dof_object (); 00532 #endif 00533 this->clear_dofs (); 00534 } 00535 00536 00537 00538 inline 00539 void DofObject::invalidate_dofs (const unsigned int sys_num) 00540 { 00541 // If the user does not specify the system number... 00542 if (sys_num >= this->n_systems()) 00543 { 00544 for (unsigned int s=0; s<this->n_systems(); s++) 00545 for (unsigned int vg=0; vg<this->n_var_groups(s); vg++) 00546 if (this->n_comp_group(s,vg)) 00547 this->set_vg_dof_base(s,vg,invalid_id); 00548 } 00549 // ...otherwise invalidate the dofs for all systems 00550 else 00551 for (unsigned int vg=0; vg<this->n_var_groups(sys_num); vg++) 00552 if (this->n_comp_group(sys_num,vg)) 00553 this->set_vg_dof_base(sys_num,vg,invalid_id); 00554 } 00555 00556 00557 00558 inline 00559 void DofObject::invalidate_id () 00560 { 00561 this->set_id (invalid_id); 00562 } 00563 00564 00565 00566 inline 00567 void DofObject::invalidate_processor_id () 00568 { 00569 this->processor_id (invalid_processor_id); 00570 } 00571 00572 00573 00574 inline 00575 void DofObject::invalidate () 00576 { 00577 this->invalidate_dofs (); 00578 this->invalidate_id (); 00579 this->invalidate_processor_id (); 00580 } 00581 00582 00583 00584 inline 00585 void DofObject::clear_dofs () 00586 { 00587 // vector swap trick to force deallocation 00588 index_buffer_t().swap(_idx_buf); 00589 00590 libmesh_assert_equal_to (this->n_systems(), 0); 00591 libmesh_assert (_idx_buf.empty()); 00592 } 00593 00594 00595 00596 inline 00597 unsigned int DofObject::n_dofs (const unsigned int s, 00598 const unsigned int var) const 00599 { 00600 libmesh_assert_less (s, this->n_systems()); 00601 00602 unsigned int num = 0; 00603 00604 // Count all variables 00605 if (var == libMesh::invalid_uint) 00606 for (unsigned int v=0; v<this->n_vars(s); v++) 00607 num += this->n_comp(s,v); 00608 00609 // Only count specified variable 00610 else 00611 num = this->n_comp(s,var); 00612 00613 return num; 00614 } 00615 00616 00617 00618 inline 00619 dof_id_type DofObject::id () const 00620 { 00621 libmesh_assert (this->valid_id()); 00622 return _id; 00623 } 00624 00625 00626 00627 inline 00628 dof_id_type & DofObject::set_id () 00629 { 00630 return _id; 00631 } 00632 00633 00634 00635 inline 00636 unique_id_type DofObject::unique_id () const 00637 { 00638 #ifdef LIBMESH_ENABLE_UNIQUE_ID 00639 libmesh_assert (this->valid_unique_id()); 00640 return _unique_id; 00641 #else 00642 return invalid_unique_id; 00643 #endif 00644 } 00645 00646 00647 00648 inline 00649 unique_id_type & DofObject::set_unique_id () 00650 { 00651 #ifdef LIBMESH_ENABLE_UNIQUE_ID 00652 return _unique_id; 00653 #else 00654 libmesh_not_implemented(); 00655 #endif 00656 } 00657 00658 00659 00660 inline 00661 bool DofObject::valid_id () const 00662 { 00663 return (DofObject::invalid_id != _id); 00664 } 00665 00666 00667 00668 inline 00669 bool DofObject::valid_unique_id () const 00670 { 00671 #ifdef LIBMESH_ENABLE_UNIQUE_ID 00672 return (DofObject::invalid_unique_id != _unique_id); 00673 #else 00674 return false; 00675 #endif 00676 } 00677 00678 00679 00680 inline 00681 processor_id_type DofObject::processor_id () const 00682 { 00683 return _processor_id; 00684 } 00685 00686 00687 00688 inline 00689 processor_id_type& DofObject::processor_id () 00690 { 00691 return _processor_id; 00692 } 00693 00694 00695 00696 inline 00697 void DofObject::processor_id (const processor_id_type pid) 00698 { 00699 this->processor_id() = pid; 00700 } 00701 00702 00703 00704 inline 00705 bool DofObject::valid_processor_id () const 00706 { 00707 return (DofObject::invalid_processor_id != _processor_id); 00708 } 00709 00710 00711 00712 inline 00713 unsigned int DofObject::n_systems () const 00714 { 00715 return _idx_buf.empty() ? 00716 0 : cast_int<unsigned int>(_idx_buf[0]); 00717 } 00718 00719 00720 00721 inline 00722 unsigned int DofObject::n_var_groups(const unsigned int s) const 00723 { 00724 libmesh_assert_less (s, this->n_systems()); 00725 00726 return (this->end_idx(s) - this->start_idx(s)) / 2; 00727 } 00728 00729 00730 00731 inline 00732 unsigned int DofObject::n_vars(const unsigned int s, 00733 const unsigned int vg) const 00734 { 00735 libmesh_assert_less (s, this->n_systems()); 00736 libmesh_assert_less (vg, this->n_var_groups(s)); 00737 00738 const unsigned int start_idx_sys = this->start_idx(s); 00739 00740 libmesh_assert_less ((start_idx_sys + 2*vg), _idx_buf.size()); 00741 00742 return (cast_int<unsigned int> 00743 (_idx_buf[start_idx_sys + 2*vg]) >> ncv_magic_exp); 00744 } 00745 00746 00747 00748 inline 00749 unsigned int DofObject::n_vars(const unsigned int s) const 00750 { 00751 libmesh_assert_less (s, this->n_systems()); 00752 00753 const unsigned int nvg = this->n_var_groups(s); 00754 00755 unsigned int val=0; 00756 00757 for (unsigned int vg=0; vg<nvg; vg++) 00758 val += this->n_vars(s,vg); 00759 00760 return val; 00761 } 00762 00763 00764 00765 00766 inline 00767 unsigned int DofObject::n_comp(const unsigned int s, 00768 const unsigned int var) const 00769 { 00770 libmesh_assert_less (s, this->n_systems()); 00771 libmesh_assert_less (var, this->n_vars(s)); 00772 00773 return this->n_comp_group(s,this->var_to_vg(s,var)); 00774 } 00775 00776 00777 00778 00779 inline 00780 unsigned int DofObject::n_comp_group(const unsigned int s, 00781 const unsigned int vg) const 00782 { 00783 libmesh_assert_less (s, this->n_systems()); 00784 libmesh_assert_less (vg, this->n_var_groups(s)); 00785 00786 const unsigned int 00787 start_idx_sys = this->start_idx(s); 00788 00789 libmesh_assert_less ((start_idx_sys + 2*vg), _idx_buf.size()); 00790 00791 return (_idx_buf[start_idx_sys + 2*vg] % ncv_magic); 00792 } 00793 00794 00795 00796 inline 00797 dof_id_type DofObject::dof_number(const unsigned int s, 00798 const unsigned int var, 00799 const unsigned int comp) const 00800 { 00801 libmesh_assert_less (s, this->n_systems()); 00802 libmesh_assert_less (var, this->n_vars(s)); 00803 libmesh_assert_less (comp, this->n_comp(s,var)); 00804 00805 const unsigned int 00806 vg = this->var_to_vg(s,var), 00807 start_idx_sys = this->start_idx(s); 00808 00809 libmesh_assert_less ((start_idx_sys + 2*vg + 1), _idx_buf.size()); 00810 00811 const dof_id_type 00812 base_idx = _idx_buf[start_idx_sys + 2*vg + 1]; 00813 00814 // if the first component is invalid, they 00815 // are all invalid 00816 if (base_idx == invalid_id) 00817 return invalid_id; 00818 00819 // otherwise the index is the first component 00820 // index augemented by the component number 00821 else 00822 { 00823 const unsigned int 00824 ncg = this->n_comp_group(s,vg), 00825 vig = this->system_var_to_vg_var(s,vg,var); 00826 00827 // std::cout << "base_idx, var, vg, vig, ncg, comp=" 00828 // << base_idx << " " 00829 // << var << " " 00830 // << vg << " " 00831 // << vig << " " 00832 // << ncg << " " 00833 // << comp << '\n'; 00834 00835 return cast_int<dof_id_type>(base_idx + vig*ncg + comp); 00836 } 00837 } 00838 00839 00840 00841 inline 00842 bool DofObject::has_dofs (const unsigned int sys) const 00843 { 00844 if (sys == libMesh::invalid_uint) 00845 { 00846 for (unsigned int s=0; s<this->n_systems(); s++) 00847 if (this->n_vars(s)) 00848 return true; 00849 } 00850 00851 else 00852 { 00853 libmesh_assert_less (sys, this->n_systems()); 00854 00855 if (this->n_vars(sys)) 00856 return true; 00857 } 00858 00859 return false; 00860 } 00861 00862 00863 00864 inline 00865 unsigned int DofObject::start_idx (const unsigned int s) const 00866 { 00867 libmesh_assert_less (s, this->n_systems()); 00868 libmesh_assert_less (s, _idx_buf.size()); 00869 00870 return cast_int<unsigned int>(_idx_buf[s]); 00871 } 00872 00873 00874 00875 inline 00876 unsigned int DofObject::end_idx (const unsigned int s) const 00877 { 00878 libmesh_assert_less (s, this->n_systems()); 00879 libmesh_assert_less (s, _idx_buf.size()); 00880 00881 return ((s+1) == this->n_systems()) ? 00882 cast_int<unsigned int>(_idx_buf.size()) : 00883 cast_int<unsigned int>(_idx_buf[s+1]); 00884 } 00885 00886 00887 00888 inline 00889 void DofObject::set_vg_dof_base(const unsigned int s, 00890 const unsigned int vg, 00891 const dof_id_type db) 00892 { 00893 libmesh_assert_less (s, this->n_systems()); 00894 libmesh_assert_less (vg, this->n_var_groups(s)); 00895 00896 const unsigned int 00897 start_idx_sys = this->start_idx(s); 00898 00899 libmesh_assert_less ((start_idx_sys + 2*vg + 1), _idx_buf.size()); 00900 00901 _idx_buf[start_idx_sys + 2*vg + 1] = db; 00902 00903 libmesh_assert_equal_to (this->vg_dof_base(s,vg), db); 00904 } 00905 00906 00907 00908 inline 00909 dof_id_type DofObject::vg_dof_base(const unsigned int s, 00910 const unsigned int vg) const 00911 { 00912 libmesh_assert_less (s, this->n_systems()); 00913 libmesh_assert_less (vg, this->n_var_groups(s)); 00914 00915 const unsigned int 00916 start_idx_sys = this->start_idx(s); 00917 00918 libmesh_assert_less ((start_idx_sys + 2*vg + 1), _idx_buf.size()); 00919 00920 // #ifdef DEBUG 00921 // std::cout << " [ "; 00922 // for (unsigned int i=0; i<_idx_buf.size(); i++) 00923 // std::cout << _idx_buf[i] << " "; 00924 // std::cout << "]\n"; 00925 // #endif 00926 00927 return _idx_buf[start_idx_sys + 2*vg + 1]; 00928 } 00929 00930 00931 00932 inline 00933 unsigned int DofObject::var_to_vg (const unsigned int s, 00934 const unsigned int var) const 00935 { 00936 const unsigned int 00937 nvg = this->n_var_groups(s); 00938 00939 for (unsigned int vg=0, vg_end=0; vg<nvg; vg++) 00940 { 00941 vg_end += this->n_vars(s,vg); 00942 if (var < vg_end) return vg; 00943 } 00944 00945 libmesh_error_msg("We'll never get here!"); 00946 return 0; 00947 } 00948 00949 00950 00951 inline 00952 unsigned int DofObject::system_var_to_vg_var (const unsigned int s, 00953 const unsigned int vg, 00954 const unsigned int var) const 00955 { 00956 unsigned int accumulated_sum=0; 00957 00958 for (unsigned int vgc=0; vgc<vg; vgc++) 00959 accumulated_sum += this->n_vars(s,vgc); 00960 00961 libmesh_assert_less_equal (accumulated_sum, var); 00962 00963 return (var - accumulated_sum); 00964 } 00965 00966 00967 } // namespace libMesh 00968 00969 00970 #endif // #ifndef LIBMESH_DOF_OBJECT_H