$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_INF_FE_H 00021 #define LIBMESH_INF_FE_H 00022 00023 #include "libmesh/libmesh_config.h" 00024 00025 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00026 00027 // Local includes 00028 #include "libmesh/fe_base.h" 00029 00030 // C++ includes 00031 #include <cstddef> 00032 00033 namespace libMesh 00034 { 00035 00036 00037 // forward declarations 00038 class Elem; 00039 class FEComputeData; 00040 00041 00042 00074 //------------------------------------------------------------- 00075 // InfFE class definition 00076 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00077 class InfFE : public FEBase 00078 { 00079 00080 /* 00081 * Protect the nested class 00082 */ 00083 protected: 00084 00095 //------------------------------------------------------------- 00096 // InfFE::Radial class definition 00097 class Radial 00098 { 00099 private: 00100 00104 Radial() {} 00105 00106 public: 00107 00112 static Real decay (const Real v); 00113 00118 static Real decay_deriv (const Real) { return -.5; } 00119 00124 static Real D (const Real v) { return (1.-v)*(1.-v)/4.; } 00125 00129 static Real D_deriv (const Real v) { return (v-1.)/2.; } 00130 00135 static Order mapping_order() { return FIRST; } 00136 00149 static unsigned int n_dofs (const Order o_radial) 00150 { return static_cast<unsigned int>(o_radial)+1; } 00151 00161 static unsigned int n_dofs_at_node (const Order o_radial, 00162 const unsigned int n_onion); 00163 00171 static unsigned int n_dofs_per_elem (const Order o_radial) 00172 { return static_cast<unsigned int>(o_radial)+1; } 00173 00174 }; 00175 00176 00177 00186 //------------------------------------------------------------- 00187 // InfFE::Base class definition 00188 class Base 00189 { 00190 private: 00191 00195 Base() {} 00196 00197 public: 00198 00204 static Elem* build_elem (const Elem* inf_elem); 00205 00211 static ElemType get_elem_type (const ElemType type); 00212 00218 static unsigned int n_base_mapping_sf (const ElemType base_elem_type, 00219 const Order base_mapping_order); 00220 00221 }; 00222 00223 00224 00225 00226 00227 00228 00229 public: 00230 00231 //------------------------------------------------------------- 00232 // InfFE continued 00233 00247 explicit 00248 InfFE(const FEType& fet); 00249 00253 ~InfFE(); 00254 00255 00256 00257 00258 00259 //------------------------------------------------------------- 00260 // The static public members for access from FEInterface etc 00272 static Real shape(const FEType& fet, 00273 const ElemType t, 00274 const unsigned int i, 00275 const Point& p); 00276 00288 static Real shape(const FEType& fet, 00289 const Elem* elem, 00290 const unsigned int i, 00291 const Point& p); 00292 00302 static void compute_data(const FEType& fe_t, 00303 const Elem* inf_elem, 00304 FEComputeData& data); 00305 00310 static unsigned int n_shape_functions (const FEType& fet, 00311 const ElemType t) 00312 { return n_dofs(fet, t); } 00313 00319 static unsigned int n_dofs(const FEType& fet, 00320 const ElemType inf_elem_type); 00321 00326 static unsigned int n_dofs_at_node(const FEType& fet, 00327 const ElemType inf_elem_type, 00328 const unsigned int n); 00329 00334 static unsigned int n_dofs_per_elem(const FEType& fet, 00335 const ElemType inf_elem_type); 00336 00340 virtual FEContinuity get_continuity() const 00341 { return C_ZERO; } // FIXME - is this true?? 00342 00347 virtual bool is_hierarchic() const 00348 { return false; } // FIXME - Inf FEs don't handle p elevation yet 00349 00357 static void nodal_soln(const FEType& fet, 00358 const Elem* elem, 00359 const std::vector<Number>& elem_soln, 00360 std::vector<Number>& nodal_soln); 00361 00377 static Point inverse_map (const Elem* elem, 00378 const Point& p, 00379 const Real tolerance = TOLERANCE, 00380 const bool secure = true, 00381 const bool interpolated = true); 00382 00383 00390 static void inverse_map (const Elem* elem, 00391 const std::vector<Point>& physical_points, 00392 std::vector<Point>& reference_points, 00393 const Real tolerance = TOLERANCE, 00394 const bool secure = true); 00395 00396 00397 //------------------------------------------------------------- 00398 // The work-horses of InfFE. These are often used during matrix assembly 00405 virtual void reinit (const Elem* elem, 00406 const std::vector<Point>* const pts = NULL, 00407 const std::vector<Real>* const weights = NULL); 00408 00414 virtual void reinit (const Elem* elem, 00415 const unsigned int side, 00416 const Real tolerance = TOLERANCE, 00417 const std::vector<Point>* const pts = NULL, 00418 const std::vector<Real>* const weights = NULL); 00419 00425 virtual void edge_reinit (const Elem* elem, 00426 const unsigned int edge, 00427 const Real tolerance = TOLERANCE, 00428 const std::vector<Point>* const pts = NULL, 00429 const std::vector<Real>* const weights = NULL); 00430 00435 virtual void side_map (const Elem* /* elem */, 00436 const Elem* /* side */, 00437 const unsigned int /* s */, 00438 const std::vector<Point>& /* reference_side_points */, 00439 std::vector<Point>& /* reference_points */) 00440 { 00441 libmesh_not_implemented(); 00442 } 00443 00455 virtual void attach_quadrature_rule (QBase* q); 00456 00461 virtual unsigned int n_shape_functions () const 00462 { return _n_total_approx_sf; } 00463 00469 virtual unsigned int n_quadrature_points () const 00470 { libmesh_assert(radial_qrule); return _n_total_qp; } 00471 00472 00473 protected: 00474 00475 //------------------------------------------------------------- 00476 // static members used by the "work-horses" 00477 00494 static Real eval(Real v, 00495 Order o_radial, 00496 unsigned int i); 00497 00503 static Real eval_deriv(Real v, 00504 Order o_radial, 00505 unsigned int i); 00506 00507 00508 00509 //------------------------------------------------------------- 00510 // Non-static members used by the "work-horses" 00515 void update_base_elem (const Elem* inf_elem); 00516 00520 virtual void init_base_shape_functions(const std::vector<Point>&, 00521 const Elem*) 00522 { libmesh_not_implemented(); } 00523 00529 void init_radial_shape_functions(const Elem* inf_elem); 00530 00537 void init_shape_functions(const Elem* inf_elem); 00538 00543 void init_face_shape_functions (const std::vector<Point>& qp, 00544 const Elem* side); 00545 00554 void combine_base_radial(const Elem* inf_elem); 00555 00567 virtual void compute_shape_functions(const Elem*, const std::vector<Point>&); 00568 00569 00570 00571 //------------------------------------------------------------- 00572 // Miscellaneous static members 00573 00578 static Point map (const Elem* inf_elem, 00579 const Point& reference_point); 00580 00587 static void compute_node_indices (const ElemType inf_elem_type, 00588 const unsigned int outer_node_index, 00589 unsigned int& base_node, 00590 unsigned int& radial_node); 00591 00600 static void compute_node_indices_fast (const ElemType inf_elem_type, 00601 const unsigned int outer_node_index, 00602 unsigned int& base_node, 00603 unsigned int& radial_node); 00604 00611 static void compute_shape_indices (const FEType& fet, 00612 const ElemType inf_elem_type, 00613 const unsigned int i, 00614 unsigned int& base_shape, 00615 unsigned int& radial_shape); 00616 00617 //-------------------------------------------------------------- 00618 // protected members, which are not to be accessed from outside 00622 std::vector<Real> dist; 00623 00630 std::vector<Real> dweightdv; 00631 00638 std::vector<Real> som; 00643 std::vector<Real> dsomdv; 00644 00649 std::vector<std::vector<Real> > mode; 00650 00655 std::vector<std::vector<Real> > dmodedv; 00656 00660 std::vector<std::vector<Real> > radial_map; 00661 00662 00666 std::vector<std::vector<Real> > dradialdv_map; 00667 00673 std::vector<Real> dphasedxi; 00674 00680 std::vector<Real> dphasedeta; 00681 00687 std::vector<Real> dphasedzeta; 00688 00689 00690 00691 00692 //-------------------------------------------------------------- 00693 // numbering scheme maps 00694 00703 std::vector<unsigned int> _radial_node_index; 00704 00713 std::vector<unsigned int> _base_node_index; 00714 00723 std::vector<unsigned int> _radial_shape_index; 00724 00733 std::vector<unsigned int> _base_shape_index; 00734 00735 00736 00737 00738 //-------------------------------------------------------------- 00739 // some more protected members 00740 00745 unsigned int _n_total_approx_sf; 00746 00751 unsigned int _n_total_qp; 00752 00757 std::vector<Real> _total_qrule_weights; 00758 00763 QBase* base_qrule; 00764 00769 QBase* radial_qrule; 00770 00775 Elem* base_elem; 00776 00783 FEBase* base_fe; 00784 00792 FEType current_fe_type; 00793 00794 00795 private: 00796 00800 virtual bool shapes_need_reinit() const; 00801 00809 static ElemType _compute_node_indices_fast_current_elem_type; 00810 00811 00812 #ifdef DEBUG 00813 00817 static bool _warned_for_nodal_soln; 00818 static bool _warned_for_shape; 00819 00820 #endif 00821 00827 template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map> 00828 friend class InfFE; 00829 00830 }; 00831 00832 00833 00834 00835 // ------------------------------------------------------------ 00836 // InfFE class inline members 00837 00838 00839 00840 00841 // ------------------------------------------------------------ 00842 // InfFE::Radial class inline members 00843 00844 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00845 inline 00846 Real InfFE<Dim,T_radial,T_map>::Radial::decay(const Real v) 00847 { 00848 switch (Dim) 00849 //TODO:[DD] What decay do i have in 2D and 1D? 00850 { 00851 case 3: 00852 return (1.-v)/2.; 00853 00854 case 2: 00855 return 0.; 00856 00857 case 1: 00858 return 0.; 00859 00860 default: 00861 libmesh_error_msg("Invalid Dim = " << Dim); 00862 } 00863 } 00864 00865 00866 00867 // ------------------------------------------------------------ 00868 // InfFE::Base class inline members 00869 00870 00871 } // namespace libMesh 00872 00873 00874 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00875 00876 00877 #endif // LIBMESH_INF_FE_H