dof_object.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 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 void set_id (const dof_id_type dofid) 00152 { this->set_id() = dofid; } 00153 00158 bool valid_id () const; 00159 00163 processor_id_type processor_id () const; 00164 00169 processor_id_type& processor_id (); 00170 00174 void processor_id (const processor_id_type pid); 00175 00180 bool valid_processor_id () const; 00181 00186 unsigned int n_systems() const; 00187 00191 void set_n_systems (const unsigned int s); 00192 00196 void add_system (); 00197 00202 unsigned int n_var_groups(const unsigned int s) const; 00203 00208 unsigned int n_vars(const unsigned int s, 00209 const unsigned int vg) const; 00210 00215 unsigned int n_vars(const unsigned int s) const; 00216 00224 void set_n_vars_per_group(const unsigned int s, 00225 const std::vector<unsigned int> &nvpg); 00226 00236 unsigned int n_comp(const unsigned int s, 00237 const unsigned int var) const; 00238 00248 unsigned int n_comp_group(const unsigned int s, 00249 const unsigned int vg) const; 00250 00255 void set_n_comp(const unsigned int s, 00256 const unsigned int var, 00257 const unsigned int ncomp); 00258 00263 void set_n_comp_group(const unsigned int s, 00264 const unsigned int vg, 00265 const unsigned int ncomp); 00266 00271 dof_id_type dof_number(const unsigned int s, 00272 const unsigned int var, 00273 const unsigned int comp) const; 00274 00279 void set_dof_number(const unsigned int s, 00280 const unsigned int var, 00281 const unsigned int comp, 00282 const dof_id_type dn); 00283 00288 bool has_dofs(const unsigned int s=libMesh::invalid_uint) const; 00289 00295 void set_vg_dof_base(const unsigned int s, 00296 const unsigned int vg, 00297 const dof_id_type db); 00298 00304 dof_id_type vg_dof_base(const unsigned int s, 00305 const unsigned int vg) const; 00306 00310 static const dof_id_type invalid_id = static_cast<dof_id_type>(-1); 00311 00316 static const processor_id_type invalid_processor_id = static_cast<processor_id_type>(-1); 00317 00322 unsigned int packed_indexing_size() const; 00323 00328 static unsigned int unpackable_indexing_size(std::vector<int>::const_iterator begin); 00329 00335 void unpack_indexing(std::vector<int>::const_iterator begin); 00336 00342 void pack_indexing(std::back_insert_iterator<std::vector<int> > target) const; 00343 00347 void debug_buffer () const; 00348 00349 private: 00350 00355 unsigned int var_to_vg (const unsigned int s, 00356 const unsigned int var) const; 00357 00362 unsigned int system_var_to_vg_var (const unsigned int s, 00363 const unsigned int vg, 00364 const unsigned int var) const; 00365 00369 dof_id_type _id; 00370 00380 processor_id_type _processor_id; 00381 00424 typedef dof_id_type index_t; 00425 typedef std::vector<index_t> index_buffer_t; 00426 index_buffer_t _idx_buf; 00427 00437 static const index_t ncv_magic = 256; // = 2^8, in case we want to manually bitshift 00438 00442 unsigned int start_idx(const unsigned int s) const; 00443 00447 unsigned int end_idx(const unsigned int s) const; 00448 00449 // methods only available for unit testing 00450 #ifdef LIBMESH_IS_UNIT_TESTING 00451 public: 00452 void set_buffer (const std::vector<dof_id_type> &buf) 00453 { _idx_buf = buf; } 00454 #endif 00455 }; 00456 00457 00458 00459 //------------------------------------------------------ 00460 // Inline functions 00461 inline 00462 DofObject::DofObject () : 00463 #ifdef LIBMESH_ENABLE_AMR 00464 old_dof_object(NULL), 00465 #endif 00466 _id (invalid_id), 00467 _processor_id (invalid_processor_id) 00468 { 00469 this->invalidate(); 00470 } 00471 00472 00473 00474 00475 00476 inline 00477 DofObject::~DofObject () 00478 { 00479 // Free all memory. 00480 #ifdef LIBMESH_ENABLE_AMR 00481 this->clear_old_dof_object (); 00482 #endif 00483 this->clear_dofs (); 00484 } 00485 00486 00487 00488 inline 00489 void DofObject::invalidate_dofs (const unsigned int sys_num) 00490 { 00491 // If the user does not specify the system number... 00492 if (sys_num >= this->n_systems()) 00493 { 00494 for (unsigned int s=0; s<this->n_systems(); s++) 00495 for (unsigned int vg=0; vg<this->n_var_groups(s); vg++) 00496 if (this->n_comp_group(s,vg)) 00497 this->set_vg_dof_base(s,vg,invalid_id); 00498 } 00499 // ...otherwise invalidate the dofs for all systems 00500 else 00501 for (unsigned int vg=0; vg<this->n_var_groups(sys_num); vg++) 00502 if (this->n_comp_group(sys_num,vg)) 00503 this->set_vg_dof_base(sys_num,vg,invalid_id); 00504 } 00505 00506 00507 00508 inline 00509 void DofObject::invalidate_id () 00510 { 00511 this->set_id (invalid_id); 00512 } 00513 00514 00515 00516 inline 00517 void DofObject::invalidate_processor_id () 00518 { 00519 this->processor_id (invalid_processor_id); 00520 } 00521 00522 00523 00524 inline 00525 void DofObject::invalidate () 00526 { 00527 this->invalidate_dofs (); 00528 this->invalidate_id (); 00529 this->invalidate_processor_id (); 00530 } 00531 00532 00533 00534 inline 00535 void DofObject::clear_dofs () 00536 { 00537 // vector swap trick to force deallocation 00538 index_buffer_t().swap(_idx_buf); 00539 00540 libmesh_assert_equal_to (this->n_systems(), 0); 00541 libmesh_assert (_idx_buf.empty()); 00542 } 00543 00544 00545 00546 inline 00547 unsigned int DofObject::n_dofs (const unsigned int s, 00548 const unsigned int var) const 00549 { 00550 libmesh_assert_less (s, this->n_systems()); 00551 00552 unsigned int num = 0; 00553 00554 // Count all variables 00555 if (var == libMesh::invalid_uint) 00556 for (unsigned int v=0; v<this->n_vars(s); v++) 00557 num += this->n_comp(s,v); 00558 00559 // Only count specified variable 00560 else 00561 num = this->n_comp(s,var); 00562 00563 return num; 00564 } 00565 00566 00567 00568 inline 00569 dof_id_type DofObject::id () const 00570 { 00571 libmesh_assert (this->valid_id()); 00572 return _id; 00573 } 00574 00575 00576 00577 inline 00578 dof_id_type & DofObject::set_id () 00579 { 00580 return _id; 00581 } 00582 00583 00584 00585 inline 00586 bool DofObject::valid_id () const 00587 { 00588 return (DofObject::invalid_id != _id); 00589 } 00590 00591 00592 inline 00593 processor_id_type DofObject::processor_id () const 00594 { 00595 return _processor_id; 00596 } 00597 00598 00599 00600 inline 00601 processor_id_type& DofObject::processor_id () 00602 { 00603 return _processor_id; 00604 } 00605 00606 00607 00608 inline 00609 void DofObject::processor_id (const processor_id_type pid) 00610 { 00611 this->processor_id() = pid; 00612 } 00613 00614 00615 00616 inline 00617 bool DofObject::valid_processor_id () const 00618 { 00619 return (DofObject::invalid_processor_id != _processor_id); 00620 } 00621 00622 00623 00624 inline 00625 unsigned int DofObject::n_systems () const 00626 { 00627 return _idx_buf.empty() ? 00628 0 : libmesh_cast_int<unsigned int>(_idx_buf[0]); 00629 } 00630 00631 00632 00633 inline 00634 unsigned int DofObject::n_var_groups(const unsigned int s) const 00635 { 00636 libmesh_assert_less (s, this->n_systems()); 00637 00638 return (this->end_idx(s) - this->start_idx(s)) / 2; 00639 } 00640 00641 00642 00643 inline 00644 unsigned int DofObject::n_vars(const unsigned int s, 00645 const unsigned int vg) const 00646 { 00647 libmesh_assert_less (s, this->n_systems()); 00648 libmesh_assert_less (vg, this->n_var_groups(s)); 00649 00650 const unsigned int start_idx_sys = this->start_idx(s); 00651 00652 libmesh_assert_less ((start_idx_sys + 2*vg), _idx_buf.size()); 00653 00654 return (libmesh_cast_int<unsigned int> 00655 (_idx_buf[start_idx_sys + 2*vg]) / ncv_magic); 00656 } 00657 00658 00659 00660 inline 00661 unsigned int DofObject::n_vars(const unsigned int s) const 00662 { 00663 libmesh_assert_less (s, this->n_systems()); 00664 00665 const unsigned int nvg = this->n_var_groups(s); 00666 00667 unsigned int val=0; 00668 00669 for (unsigned int vg=0; vg<nvg; vg++) 00670 val += this->n_vars(s,vg); 00671 00672 return val; 00673 } 00674 00675 00676 00677 00678 inline 00679 unsigned int DofObject::n_comp(const unsigned int s, 00680 const unsigned int var) const 00681 { 00682 libmesh_assert_less (s, this->n_systems()); 00683 libmesh_assert_less (var, this->n_vars(s)); 00684 00685 return this->n_comp_group(s,this->var_to_vg(s,var)); 00686 } 00687 00688 00689 00690 00691 inline 00692 unsigned int DofObject::n_comp_group(const unsigned int s, 00693 const unsigned int vg) const 00694 { 00695 libmesh_assert_less (s, this->n_systems()); 00696 libmesh_assert_less (vg, this->n_var_groups(s)); 00697 00698 const unsigned int 00699 start_idx_sys = this->start_idx(s); 00700 00701 libmesh_assert_less ((start_idx_sys + 2*vg), _idx_buf.size()); 00702 00703 return (_idx_buf[start_idx_sys + 2*vg] % ncv_magic); 00704 } 00705 00706 00707 00708 inline 00709 dof_id_type DofObject::dof_number(const unsigned int s, 00710 const unsigned int var, 00711 const unsigned int comp) const 00712 { 00713 libmesh_assert_less (s, this->n_systems()); 00714 libmesh_assert_less (var, this->n_vars(s)); 00715 libmesh_assert_less (comp, this->n_comp(s,var)); 00716 00717 const unsigned int 00718 vg = this->var_to_vg(s,var), 00719 start_idx_sys = this->start_idx(s); 00720 00721 libmesh_assert_less ((start_idx_sys + 2*vg + 1), _idx_buf.size()); 00722 00723 const dof_id_type 00724 base_idx = _idx_buf[start_idx_sys + 2*vg + 1]; 00725 00726 // if the first component is invalid, they 00727 // are all invalid 00728 if (base_idx == invalid_id) 00729 return invalid_id; 00730 00731 // otherwise the index is the first component 00732 // index augemented by the component number 00733 else 00734 { 00735 const unsigned int 00736 ncg = this->n_comp_group(s,vg), 00737 vig = this->system_var_to_vg_var(s,vg,var); 00738 00739 // std::cout << "base_idx, var, vg, vig, ncg, comp=" 00740 // << base_idx << " " 00741 // << var << " " 00742 // << vg << " " 00743 // << vig << " " 00744 // << ncg << " " 00745 // << comp << '\n'; 00746 00747 return libmesh_cast_int<dof_id_type>(base_idx + vig*ncg + comp); 00748 } 00749 } 00750 00751 00752 00753 inline 00754 bool DofObject::has_dofs (const unsigned int sys) const 00755 { 00756 if (sys == libMesh::invalid_uint) 00757 { 00758 for (unsigned int s=0; s<this->n_systems(); s++) 00759 if (this->n_vars(s)) 00760 return true; 00761 } 00762 00763 else 00764 { 00765 libmesh_assert_less (sys, this->n_systems()); 00766 00767 if (this->n_vars(sys)) 00768 return true; 00769 } 00770 00771 return false; 00772 } 00773 00774 00775 00776 inline 00777 unsigned int DofObject::start_idx (const unsigned int s) const 00778 { 00779 libmesh_assert_less (s, this->n_systems()); 00780 libmesh_assert_less (s, _idx_buf.size()); 00781 00782 return libmesh_cast_int<unsigned int>(_idx_buf[s]); 00783 } 00784 00785 00786 00787 inline 00788 unsigned int DofObject::end_idx (const unsigned int s) const 00789 { 00790 libmesh_assert_less (s, this->n_systems()); 00791 libmesh_assert_less (s, _idx_buf.size()); 00792 00793 return ((s+1) == this->n_systems()) ? 00794 libmesh_cast_int<unsigned int>(_idx_buf.size()) : 00795 libmesh_cast_int<unsigned int>(_idx_buf[s+1]); 00796 } 00797 00798 00799 00800 inline 00801 void DofObject::set_vg_dof_base(const unsigned int s, 00802 const unsigned int vg, 00803 const dof_id_type db) 00804 { 00805 libmesh_assert_less (s, this->n_systems()); 00806 libmesh_assert_less (vg, this->n_var_groups(s)); 00807 00808 const unsigned int 00809 start_idx_sys = this->start_idx(s); 00810 00811 libmesh_assert_less ((start_idx_sys + 2*vg + 1), _idx_buf.size()); 00812 00813 _idx_buf[start_idx_sys + 2*vg + 1] = db; 00814 00815 libmesh_assert_equal_to (this->vg_dof_base(s,vg), db); 00816 } 00817 00818 00819 00820 inline 00821 dof_id_type DofObject::vg_dof_base(const unsigned int s, 00822 const unsigned int vg) const 00823 { 00824 libmesh_assert_less (s, this->n_systems()); 00825 libmesh_assert_less (vg, this->n_var_groups(s)); 00826 00827 const unsigned int 00828 start_idx_sys = this->start_idx(s); 00829 00830 libmesh_assert_less ((start_idx_sys + 2*vg + 1), _idx_buf.size()); 00831 00832 // #ifdef DEBUG 00833 // std::cout << " [ "; 00834 // for (unsigned int i=0; i<_idx_buf.size(); i++) 00835 // std::cout << _idx_buf[i] << " "; 00836 // std::cout << "]\n"; 00837 // #endif 00838 00839 return _idx_buf[start_idx_sys + 2*vg + 1]; 00840 } 00841 00842 00843 00844 inline 00845 unsigned int DofObject::var_to_vg (const unsigned int s, 00846 const unsigned int var) const 00847 { 00848 const unsigned int 00849 nvg = this->n_var_groups(s); 00850 00851 for (unsigned int vg=0, vg_end=0; vg<nvg; vg++) 00852 { 00853 vg_end += this->n_vars(s,vg); 00854 if (var < vg_end) return vg; 00855 } 00856 00857 // we should never get here 00858 libmesh_error(); 00859 return 0; 00860 } 00861 00862 00863 00864 inline 00865 unsigned int DofObject::system_var_to_vg_var (const unsigned int s, 00866 const unsigned int vg, 00867 const unsigned int var) const 00868 { 00869 unsigned int accumulated_sum=0; 00870 00871 for (unsigned int vgc=0; vgc<vg; vgc++) 00872 accumulated_sum += this->n_vars(s,vgc); 00873 00874 libmesh_assert_less_equal (accumulated_sum, var); 00875 00876 return (var - accumulated_sum); 00877 } 00878 00879 00880 } // namespace libMesh 00881 00882 00883 #endif // #ifndef LIBMESH_DOF_OBJECT_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: