inf_fe_boundary.C
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 // C++ includes 00021 00022 00023 // Local includes 00024 #include "libmesh/libmesh_config.h" 00025 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00026 #include "libmesh/inf_fe.h" 00027 #include "libmesh/inf_fe_macro.h" 00028 #include "libmesh/quadrature.h" 00029 #include "libmesh/elem.h" 00030 00031 namespace libMesh 00032 { 00033 00034 00035 00036 00037 //------------------------------------------------------- 00038 // Method for 2D, 3D -- see inf_fe_1D.C for a 1D version of this 00039 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base> 00040 void InfFE<Dim,T_radial,T_base>::reinit(const Elem* inf_elem, 00041 const unsigned int s, 00042 const Real tolerance, 00043 const std::vector<Point>* const pts, 00044 const std::vector<Real>* const weights) 00045 { 00046 if (weights != NULL) 00047 { 00048 libMesh::err << "ERROR: User-specified weights for infinite elements " 00049 << "are not implemented!" << std::endl; 00050 libmesh_not_implemented(); 00051 } 00052 00053 if (pts != NULL) 00054 { 00055 libMesh::err << "ERROR: User-specified points for infinite elements " 00056 << "are not implemented!" << std::endl; 00057 libmesh_not_implemented(); 00058 } 00059 00060 // We don't do this for 1D elements! 00061 libmesh_assert_not_equal_to (Dim, 1); 00062 00063 libmesh_assert(inf_elem); 00064 libmesh_assert(qrule); 00065 00066 // Don't do this for the base 00067 libmesh_assert_not_equal_to (s, 0); 00068 00069 // Build the side of interest 00070 const AutoPtr<Elem> side(inf_elem->build_side(s)); 00071 00072 // set the element type 00073 elem_type = inf_elem->type(); 00074 00075 // eventually initialize radial quadrature rule 00076 bool radial_qrule_initialized = false; 00077 00078 if (current_fe_type.radial_order != fe_type.radial_order) 00079 { 00080 current_fe_type.radial_order = fe_type.radial_order; 00081 radial_qrule->init(EDGE2, inf_elem->p_level()); 00082 radial_qrule_initialized = true; 00083 } 00084 00085 // Initialize the face shape functions 00086 if (this->get_type() != inf_elem->type() || 00087 base_fe->shapes_need_reinit() || 00088 radial_qrule_initialized) 00089 this->init_face_shape_functions (qrule->get_points(), side.get()); 00090 00091 00092 // compute the face map 00093 this->_fe_map->compute_face_map(this->dim, _total_qrule_weights, side.get()); 00094 00095 // make a copy of the Jacobian for integration 00096 const std::vector<Real> JxW_int(this->_fe_map->get_JxW()); 00097 00098 // Find where the integration points are located on the 00099 // full element. 00100 std::vector<Point> qp; this->inverse_map (inf_elem, this->_fe_map->get_xyz(), 00101 qp, tolerance); 00102 00103 // compute the shape function and derivative values 00104 // at the points qp 00105 this->reinit (inf_elem, &qp); 00106 00107 // copy back old data 00108 this->_fe_map->get_JxW() = JxW_int; 00109 00110 } 00111 00112 00113 00114 //------------------------------------------------------- 00115 // Method for 2D, 3D -- see inf_fe_1D.C for a 1D version of this 00116 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base> 00117 void InfFE<Dim,T_radial,T_base>::edge_reinit(const Elem*, 00118 const unsigned int, 00119 const Real, 00120 const std::vector<Point>* const pts, 00121 const std::vector<Real>* const /*weights*/) 00122 { 00123 // We don't do this for 1D elements! 00124 //libmesh_assert_not_equal_to (Dim, 1); 00125 00126 libMesh::err << "ERROR: Edge conditions for infinite elements " 00127 << "not implemented!" << std::endl; 00128 libmesh_error(); 00129 00130 if (pts != NULL) 00131 { 00132 libMesh::err << "ERROR: User-specified points for infinite elements " 00133 << "not implemented!" << std::endl; 00134 libmesh_error(); 00135 } 00136 } 00137 00138 00139 00140 00141 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base> 00142 void InfFE<Dim,T_radial,T_base>::init_face_shape_functions(const std::vector<Point>&, 00143 const Elem* inf_side) 00144 { 00145 libmesh_assert(inf_side); 00146 00147 // Currently, this makes only sense in 3-D! 00148 libmesh_assert_equal_to (Dim, 3); 00149 00150 // Initialiize the radial shape functions 00151 this->init_radial_shape_functions(inf_side); 00152 00153 // Initialize the base shape functions 00154 this->update_base_elem(inf_side); 00155 00156 // Initialize the base quadratur rule 00157 base_qrule->init(base_elem->type(), inf_side->p_level()); 00158 00159 // base_fe still corresponds to the (dim-1)-dimensional base of the InfFE object, 00160 // so update the fe_base. 00161 { 00162 libmesh_assert_equal_to (Dim, 3); 00163 00164 AutoPtr<FEBase> ap_fb(FEBase::build(Dim-2, this->fe_type)); 00165 if (base_fe != NULL) 00166 delete base_fe; 00167 base_fe = ap_fb.release(); 00168 base_fe->attach_quadrature_rule(base_qrule); 00169 } 00170 00171 // initialize the shape functions on the base 00172 base_fe->init_base_shape_functions(base_fe->qrule->get_points(), 00173 base_elem); 00174 00175 // the number of quadrature points 00176 const unsigned int n_radial_qp = 00177 libmesh_cast_int<unsigned int>(som.size()); 00178 const unsigned int n_base_qp = base_qrule->n_points(); 00179 const unsigned int n_total_qp = n_radial_qp * n_base_qp; 00180 00181 // the quadratur weigths 00182 _total_qrule_weights.resize(n_total_qp); 00183 00184 // now inite the shapes for boundary work 00185 { 00186 00187 // The element type and order to use in the base map 00188 const Order base_mapping_order ( base_elem->default_order() ); 00189 const ElemType base_mapping_elem_type ( base_elem->type() ); 00190 00191 // the number of mapping shape functions 00192 // (Lagrange shape functions are used for mapping in the base) 00193 const unsigned int n_radial_mapping_sf = 00194 libmesh_cast_int<unsigned int>(radial_map.size()); 00195 const unsigned int n_base_mapping_shape_functions = Base::n_base_mapping_sf(base_mapping_elem_type, 00196 base_mapping_order); 00197 00198 const unsigned int n_total_mapping_shape_functions = 00199 n_radial_mapping_sf * n_base_mapping_shape_functions; 00200 00201 00202 // initialize the node and shape numbering maps 00203 { 00204 _radial_node_index.resize (n_total_mapping_shape_functions); 00205 _base_node_index.resize (n_total_mapping_shape_functions); 00206 00207 const ElemType inf_face_elem_type (inf_side->type()); 00208 00209 // fill the node index map 00210 for (unsigned int n=0; n<n_total_mapping_shape_functions; n++) 00211 { 00212 compute_node_indices (inf_face_elem_type, 00213 n, 00214 _base_node_index[n], 00215 _radial_node_index[n]); 00216 00217 libmesh_assert_less (_base_node_index[n], n_base_mapping_shape_functions); 00218 libmesh_assert_less (_radial_node_index[n], n_radial_mapping_sf); 00219 } 00220 00221 } 00222 00223 // rezise map data fields 00224 { 00225 std::vector<std::vector<Real> >& psi_map = this->_fe_map->get_psi(); 00226 std::vector<std::vector<Real> >& dpsidxi_map = this->_fe_map->get_dpsidxi(); 00227 std::vector<std::vector<Real> >& d2psidxi2_map = this->_fe_map->get_d2psidxi2(); 00228 psi_map.resize (n_total_mapping_shape_functions); 00229 dpsidxi_map.resize (n_total_mapping_shape_functions); 00230 d2psidxi2_map.resize (n_total_mapping_shape_functions); 00231 00232 // if (Dim == 3) 00233 { 00234 std::vector<std::vector<Real> >& dpsideta_map = this->_fe_map->get_dpsideta(); 00235 std::vector<std::vector<Real> >& d2psidxideta_map = this->_fe_map->get_d2psidxideta(); 00236 std::vector<std::vector<Real> >& d2psideta2_map = this->_fe_map->get_d2psideta2(); 00237 dpsideta_map.resize (n_total_mapping_shape_functions); 00238 d2psidxideta_map.resize (n_total_mapping_shape_functions); 00239 d2psideta2_map.resize (n_total_mapping_shape_functions); 00240 } 00241 00242 for (unsigned int i=0; i<n_total_mapping_shape_functions; i++) 00243 { 00244 psi_map[i].resize (n_total_qp); 00245 dpsidxi_map[i].resize (n_total_qp); 00246 d2psidxi2_map[i].resize (n_total_qp); 00247 00248 // if (Dim == 3) 00249 { 00250 std::vector<std::vector<Real> >& dpsideta_map = this->_fe_map->get_dpsideta(); 00251 std::vector<std::vector<Real> >& d2psidxideta_map = this->_fe_map->get_d2psidxideta(); 00252 std::vector<std::vector<Real> >& d2psideta2_map = this->_fe_map->get_d2psideta2(); 00253 dpsideta_map[i].resize (n_total_qp); 00254 d2psidxideta_map[i].resize (n_total_qp); 00255 d2psideta2_map[i].resize (n_total_qp); 00256 } 00257 } 00258 } 00259 00260 00261 // compute shape maps 00262 { 00263 const std::vector<std::vector<Real> >& S_map = (base_fe->get_fe_map()).get_phi_map(); 00264 const std::vector<std::vector<Real> >& Ss_map = (base_fe->get_fe_map()).get_dphidxi_map(); 00265 00266 std::vector<std::vector<Real> >& psi_map = this->_fe_map->get_psi(); 00267 std::vector<std::vector<Real> >& dpsidxi_map = this->_fe_map->get_dpsidxi(); 00268 std::vector<std::vector<Real> >& dpsideta_map = this->_fe_map->get_dpsideta(); 00269 00270 for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qp's 00271 for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qp's 00272 for (unsigned int ti=0; ti<n_total_mapping_shape_functions; ti++) // over all mapping shapes 00273 { 00274 // let the index vectors take care of selecting the appropriate base/radial mapping shape 00275 const unsigned int bi = _base_node_index [ti]; 00276 const unsigned int ri = _radial_node_index[ti]; 00277 psi_map [ti][bp+rp*n_base_qp] = S_map [bi][bp] * radial_map [ri][rp]; 00278 dpsidxi_map [ti][bp+rp*n_base_qp] = Ss_map[bi][bp] * radial_map [ri][rp]; 00279 dpsideta_map [ti][bp+rp*n_base_qp] = S_map [bi][bp] * dradialdv_map[ri][rp]; 00280 00281 // second derivatives are not implemented for infinite elements 00282 // d2psidxi2_map [ti][bp+rp*n_base_qp] = 0.; 00283 // d2psidxideta_map [ti][bp+rp*n_base_qp] = 0.; 00284 // d2psideta2_map [ti][bp+rp*n_base_qp] = 0.; 00285 } 00286 00287 } 00288 00289 } 00290 00291 // quadrature rule weights 00292 { 00293 const std::vector<Real>& radial_qw = radial_qrule->get_weights(); 00294 const std::vector<Real>& base_qw = base_qrule->get_weights(); 00295 00296 libmesh_assert_equal_to (radial_qw.size(), n_radial_qp); 00297 libmesh_assert_equal_to (base_qw.size(), n_base_qp); 00298 00299 for (unsigned int rp=0; rp<n_radial_qp; rp++) 00300 for (unsigned int bp=0; bp<n_base_qp; bp++) 00301 { 00302 _total_qrule_weights[ bp+rp*n_base_qp ] = radial_qw[rp] * base_qw[bp]; 00303 } 00304 } 00305 00306 } 00307 00308 00309 00310 00311 //-------------------------------------------------------------- 00312 // Explicit instantiations - doesn't make sense in 1D, but as 00313 // long as we only return errors, we are fine... ;-) 00314 //#include "libmesh/inf_fe_instantiate_1D.h" 00315 //#include "libmesh/inf_fe_instantiate_2D.h" 00316 //#include "libmesh/inf_fe_instantiate_3D.h" 00317 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,reinit(const Elem*,const unsigned int, const Real, const std::vector<Point>* const, const std::vector<Real>* const)); 00318 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,reinit(const Elem*,const unsigned int, const Real, const std::vector<Point>* const, const std::vector<Real>* const)); 00319 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,reinit(const Elem*,const unsigned int, const Real, const std::vector<Point>* const, const std::vector<Real>* const)); 00320 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,edge_reinit(const Elem*,const unsigned int, const Real, const std::vector<Point>* const, const std::vector<Real>* const)); 00321 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,edge_reinit(const Elem*,const unsigned int, const Real, const std::vector<Point>* const, const std::vector<Real>* const)); 00322 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,edge_reinit(const Elem*,const unsigned int, const Real, const std::vector<Point>* const, const std::vector<Real>* const)); 00323 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,init_face_shape_functions(const std::vector<Point>&,const Elem*)); 00324 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,init_face_shape_functions(const std::vector<Point>&,const Elem*)); 00325 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,init_face_shape_functions(const std::vector<Point>&,const Elem*)); 00326 00327 } // namespace libMesh 00328 00329 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00330
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: