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 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00022 #include <cmath> // for std::sqrt 00023 00024 00025 // Local includes 00026 #include "libmesh/libmesh_common.h" 00027 #include "libmesh/fe.h" 00028 #include "libmesh/quadrature.h" 00029 #include "libmesh/elem.h" 00030 #include "libmesh/libmesh_logging.h" 00031 00032 namespace libMesh 00033 { 00034 00035 //------------------------------------------------------- 00036 // Full specializations for useless methods in 0D, 1D 00037 #define REINIT_ERROR(_dim, _type, _func) \ 00038 template <> \ 00039 void FE<_dim,_type>::_func(const Elem*, \ 00040 const unsigned int, \ 00041 const Real, \ 00042 const std::vector<Point>* const, \ 00043 const std::vector<Real>* const) \ 00044 { \ 00045 libMesh::err << "ERROR: This method makes no sense for low-D elements!" \ 00046 << std::endl; \ 00047 libmesh_error(); \ 00048 } 00049 00050 #define SIDEMAP_ERROR(_dim, _type, _func) \ 00051 template <> \ 00052 void FE<_dim,_type>::_func(const Elem*, \ 00053 const Elem*, \ 00054 const unsigned int, \ 00055 const std::vector<Point>&, \ 00056 std::vector<Point>&) \ 00057 { \ 00058 libMesh::err << "ERROR: This method makes no sense for low-D elements!" \ 00059 << std::endl; \ 00060 libmesh_error(); \ 00061 } 00062 00063 #define FACE_EDGE_SHAPE_ERROR(_dim, _func) \ 00064 template <> \ 00065 void FEMap::_func<_dim>(const std::vector<Point>&, \ 00066 const Elem* ) \ 00067 { \ 00068 libMesh::err << "ERROR: This method makes no sense for low-D elements!" \ 00069 << std::endl; \ 00070 libmesh_error(); \ 00071 } 00072 00073 00074 REINIT_ERROR(0, CLOUGH, reinit) 00075 REINIT_ERROR(0, CLOUGH, edge_reinit) 00076 SIDEMAP_ERROR(0, CLOUGH, side_map) 00077 REINIT_ERROR(0, HERMITE, reinit) 00078 REINIT_ERROR(0, HERMITE, edge_reinit) 00079 SIDEMAP_ERROR(0, HERMITE, side_map) 00080 REINIT_ERROR(0, HIERARCHIC, reinit) 00081 REINIT_ERROR(0, HIERARCHIC, edge_reinit) 00082 SIDEMAP_ERROR(0, HIERARCHIC, side_map) 00083 REINIT_ERROR(0, L2_HIERARCHIC, reinit) 00084 REINIT_ERROR(0, L2_HIERARCHIC, edge_reinit) 00085 SIDEMAP_ERROR(0, L2_HIERARCHIC, side_map) 00086 REINIT_ERROR(0, LAGRANGE, reinit) 00087 REINIT_ERROR(0, LAGRANGE, edge_reinit) 00088 SIDEMAP_ERROR(0, LAGRANGE, side_map) 00089 REINIT_ERROR(0, LAGRANGE_VEC, reinit) 00090 REINIT_ERROR(0, LAGRANGE_VEC, edge_reinit) 00091 SIDEMAP_ERROR(0, LAGRANGE_VEC, side_map) 00092 REINIT_ERROR(0, L2_LAGRANGE, reinit) 00093 REINIT_ERROR(0, L2_LAGRANGE, edge_reinit) 00094 SIDEMAP_ERROR(0, L2_LAGRANGE, side_map) 00095 REINIT_ERROR(0, MONOMIAL, reinit) 00096 REINIT_ERROR(0, MONOMIAL, edge_reinit) 00097 SIDEMAP_ERROR(0, MONOMIAL, side_map) 00098 REINIT_ERROR(0, SCALAR, reinit) 00099 REINIT_ERROR(0, SCALAR, edge_reinit) 00100 SIDEMAP_ERROR(0, SCALAR, side_map) 00101 REINIT_ERROR(0, XYZ, reinit) 00102 REINIT_ERROR(0, XYZ, edge_reinit) 00103 SIDEMAP_ERROR(0, XYZ, side_map) 00104 REINIT_ERROR(0, NEDELEC_ONE, reinit) 00105 REINIT_ERROR(0, NEDELEC_ONE, edge_reinit) 00106 SIDEMAP_ERROR(0, NEDELEC_ONE, side_map) 00107 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00108 REINIT_ERROR(0, BERNSTEIN, reinit) 00109 REINIT_ERROR(0, BERNSTEIN, edge_reinit) 00110 SIDEMAP_ERROR(0, BERNSTEIN, side_map) 00111 REINIT_ERROR(0, SZABAB, reinit) 00112 REINIT_ERROR(0, SZABAB, edge_reinit) 00113 SIDEMAP_ERROR(0, SZABAB, side_map) 00114 #endif 00115 00116 REINIT_ERROR(1, CLOUGH, edge_reinit) 00117 REINIT_ERROR(1, HERMITE, edge_reinit) 00118 REINIT_ERROR(1, HIERARCHIC, edge_reinit) 00119 REINIT_ERROR(1, L2_HIERARCHIC, edge_reinit) 00120 REINIT_ERROR(1, LAGRANGE, edge_reinit) 00121 REINIT_ERROR(1, LAGRANGE_VEC, edge_reinit) 00122 REINIT_ERROR(1, L2_LAGRANGE, edge_reinit) 00123 REINIT_ERROR(1, XYZ, edge_reinit) 00124 REINIT_ERROR(1, MONOMIAL, edge_reinit) 00125 REINIT_ERROR(1, SCALAR, edge_reinit) 00126 REINIT_ERROR(1, NEDELEC_ONE, reinit) 00127 REINIT_ERROR(1, NEDELEC_ONE, edge_reinit) 00128 SIDEMAP_ERROR(1, NEDELEC_ONE, side_map) 00129 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00130 REINIT_ERROR(1, BERNSTEIN, edge_reinit) 00131 REINIT_ERROR(1, SZABAB, edge_reinit) 00132 #endif 00133 00134 00135 //------------------------------------------------------- 00136 // Methods for 2D, 3D 00137 template <unsigned int Dim, FEFamily T> 00138 void FE<Dim,T>::reinit(const Elem* elem, 00139 const unsigned int s, 00140 const Real /* tolerance */, 00141 const std::vector<Point>* const pts, 00142 const std::vector<Real>* const weights) 00143 { 00144 libmesh_assert(elem); 00145 libmesh_assert (this->qrule != NULL || pts != NULL); 00146 // We now do this for 1D elements! 00147 // libmesh_assert_not_equal_to (Dim, 1); 00148 00149 // Build the side of interest 00150 const AutoPtr<Elem> side(elem->build_side(s)); 00151 00152 // Find the max p_level to select 00153 // the right quadrature rule for side integration 00154 unsigned int side_p_level = elem->p_level(); 00155 if (elem->neighbor(s) != NULL) 00156 side_p_level = std::max(side_p_level, elem->neighbor(s)->p_level()); 00157 00158 // Initialize the shape functions at the user-specified 00159 // points 00160 if (pts != NULL) 00161 { 00162 // The shape functions do not correspond to the qrule 00163 this->shapes_on_quadrature = false; 00164 00165 // Initialize the face shape functions 00166 this->_fe_map->template init_face_shape_functions<Dim>(*pts, side.get()); 00167 00168 // Compute the Jacobian*Weight on the face for integration 00169 if (weights != NULL) 00170 { 00171 this->_fe_map->compute_face_map (Dim, *weights, side.get()); 00172 } 00173 else 00174 { 00175 std::vector<Real> dummy_weights (pts->size(), 1.); 00176 this->_fe_map->compute_face_map (Dim, dummy_weights, side.get()); 00177 } 00178 } 00179 // If there are no user specified points, we use the 00180 // quadrature rule 00181 else 00182 { 00183 // initialize quadrature rule 00184 this->qrule->init(side->type(), side_p_level); 00185 00186 if(this->qrule->shapes_need_reinit()) 00187 this->shapes_on_quadrature = false; 00188 00189 // FIXME - could this break if the same FE object was used 00190 // for both volume and face integrals? - RHS 00191 // We might not need to reinitialize the shape functions 00192 if ((this->get_type() != elem->type()) || 00193 (side->type() != last_side) || 00194 (this->get_p_level() != side_p_level) || 00195 this->shapes_need_reinit() || 00196 !this->shapes_on_quadrature) 00197 { 00198 // Set the element type and p_level 00199 this->elem_type = elem->type(); 00200 00201 // Set the last_side 00202 last_side = side->type(); 00203 00204 // Set the last p level 00205 this->_p_level = side_p_level; 00206 00207 // Initialize the face shape functions 00208 this->_fe_map->template init_face_shape_functions<Dim>(this->qrule->get_points(), side.get()); 00209 } 00210 00211 // Compute the Jacobian*Weight on the face for integration 00212 this->_fe_map->compute_face_map (Dim, this->qrule->get_weights(), side.get()); 00213 00214 // The shape functions correspond to the qrule 00215 this->shapes_on_quadrature = true; 00216 } 00217 00218 // make a copy of the Jacobian for integration 00219 const std::vector<Real> JxW_int(this->_fe_map->get_JxW()); 00220 00221 // make a copy of shape on quadrature info 00222 bool shapes_on_quadrature_side = this->shapes_on_quadrature; 00223 00224 // Find where the integration points are located on the 00225 // full element. 00226 const std::vector<Point>* ref_qp; 00227 if (pts != NULL) 00228 ref_qp = pts; 00229 else 00230 ref_qp = &this->qrule->get_points(); 00231 00232 std::vector<Point> qp; 00233 this->side_map(elem, side.get(), s, *ref_qp, qp); 00234 00235 // compute the shape function and derivative values 00236 // at the points qp 00237 this->reinit (elem, &qp); 00238 00239 this->shapes_on_quadrature = shapes_on_quadrature_side; 00240 00241 // copy back old data 00242 this->_fe_map->get_JxW() = JxW_int; 00243 } 00244 00245 00246 00247 template <unsigned int Dim, FEFamily T> 00248 void FE<Dim,T>::edge_reinit(const Elem* elem, 00249 const unsigned int e, 00250 const Real tolerance, 00251 const std::vector<Point>* const pts, 00252 const std::vector<Real>* const weights) 00253 { 00254 libmesh_assert(elem); 00255 libmesh_assert (this->qrule != NULL || pts != NULL); 00256 // We don't do this for 1D elements! 00257 libmesh_assert_not_equal_to (Dim, 1); 00258 00259 // Build the side of interest 00260 const AutoPtr<Elem> edge(elem->build_edge(e)); 00261 00262 // Initialize the shape functions at the user-specified 00263 // points 00264 if (pts != NULL) 00265 { 00266 // The shape functions do not correspond to the qrule 00267 this->shapes_on_quadrature = false; 00268 00269 // Initialize the edge shape functions 00270 this->_fe_map->template init_edge_shape_functions<Dim> (*pts, edge.get()); 00271 00272 // Compute the Jacobian*Weight on the face for integration 00273 if (weights != NULL) 00274 { 00275 this->_fe_map->compute_edge_map (Dim, *weights, edge.get()); 00276 } 00277 else 00278 { 00279 std::vector<Real> dummy_weights (pts->size(), 1.); 00280 this->_fe_map->compute_edge_map (Dim, dummy_weights, edge.get()); 00281 } 00282 } 00283 // If there are no user specified points, we use the 00284 // quadrature rule 00285 else 00286 { 00287 // initialize quadrature rule 00288 this->qrule->init(edge->type(), elem->p_level()); 00289 00290 if(this->qrule->shapes_need_reinit()) 00291 this->shapes_on_quadrature = false; 00292 00293 // We might not need to reinitialize the shape functions 00294 if ((this->get_type() != elem->type()) || 00295 (edge->type() != static_cast<int>(last_edge)) || // Comparison between enum and unsigned, cast the unsigned to int 00296 this->shapes_need_reinit() || 00297 !this->shapes_on_quadrature) 00298 { 00299 // Set the element type 00300 this->elem_type = elem->type(); 00301 00302 // Set the last_edge 00303 last_edge = edge->type(); 00304 00305 // Initialize the edge shape functions 00306 this->_fe_map->template init_edge_shape_functions<Dim> (this->qrule->get_points(), edge.get()); 00307 } 00308 00309 // Compute the Jacobian*Weight on the face for integration 00310 this->_fe_map->compute_edge_map (Dim, this->qrule->get_weights(), edge.get()); 00311 00312 // The shape functions correspond to the qrule 00313 this->shapes_on_quadrature = true; 00314 } 00315 00316 // make a copy of the Jacobian for integration 00317 const std::vector<Real> JxW_int(this->_fe_map->get_JxW()); 00318 00319 // Find where the integration points are located on the 00320 // full element. 00321 std::vector<Point> qp; 00322 this->inverse_map (elem, this->_fe_map->get_xyz(), qp, tolerance); 00323 00324 // compute the shape function and derivative values 00325 // at the points qp 00326 this->reinit (elem, &qp); 00327 00328 // copy back old data 00329 this->_fe_map->get_JxW() = JxW_int; 00330 } 00331 00332 template <unsigned int Dim, FEFamily T> 00333 void FE<Dim,T>::side_map (const Elem* elem, 00334 const Elem* side, 00335 const unsigned int s, 00336 const std::vector<Point>& reference_side_points, 00337 std::vector<Point>& reference_points) 00338 { 00339 unsigned int side_p_level = elem->p_level(); 00340 if (elem->neighbor(s) != NULL) 00341 side_p_level = std::max(side_p_level, elem->neighbor(s)->p_level()); 00342 00343 if (side->type() != last_side || 00344 side_p_level != this->_p_level || 00345 !this->shapes_on_quadrature) 00346 { 00347 // Set the element type 00348 this->elem_type = elem->type(); 00349 this->_p_level = side_p_level; 00350 00351 // Set the last_side 00352 last_side = side->type(); 00353 00354 // Initialize the face shape functions 00355 this->_fe_map->template init_face_shape_functions<Dim>(reference_side_points, side); 00356 } 00357 00358 const unsigned int n_points = 00359 libmesh_cast_int<unsigned int>(reference_side_points.size()); 00360 reference_points.resize(n_points); 00361 for (unsigned int i = 0; i < n_points; i++) 00362 reference_points[i].zero(); 00363 00364 std::vector<unsigned int> elem_nodes_map; 00365 elem_nodes_map.resize(side->n_nodes()); 00366 for (unsigned int j = 0; j < side->n_nodes(); j++) 00367 for (unsigned int i = 0; i < elem->n_nodes(); i++) 00368 if (side->node(j) == elem->node(i)) 00369 elem_nodes_map[j] = i; 00370 std::vector<Point> refspace_nodes; 00371 this->get_refspace_nodes(elem->type(), refspace_nodes); 00372 00373 const std::vector<std::vector<Real> >& psi_map = this->_fe_map->get_psi(); 00374 00375 for (unsigned int i=0; i<psi_map.size(); i++) // sum over the nodes 00376 { 00377 const Point& side_node = refspace_nodes[elem_nodes_map[i]]; 00378 for (unsigned int p=0; p<n_points; p++) 00379 reference_points[p].add_scaled (side_node, psi_map[i][p]); 00380 } 00381 } 00382 00383 template<unsigned int Dim> 00384 void FEMap::init_face_shape_functions(const std::vector<Point>& qp, 00385 const Elem* side) 00386 { 00387 libmesh_assert(side); 00388 00392 START_LOG("init_face_shape_functions()", "FEMap"); 00393 00394 // The element type and order to use in 00395 // the map 00396 const Order mapping_order (side->default_order()); 00397 const ElemType mapping_elem_type (side->type()); 00398 00399 // The number of quadrature points. 00400 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qp.size()); 00401 00402 const unsigned int n_mapping_shape_functions = 00403 FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type, 00404 mapping_order); 00405 00406 // resize the vectors to hold current data 00407 // Psi are the shape functions used for the FE mapping 00408 this->psi_map.resize (n_mapping_shape_functions); 00409 00410 if (Dim > 1) 00411 { 00412 this->dpsidxi_map.resize (n_mapping_shape_functions); 00413 this->d2psidxi2_map.resize (n_mapping_shape_functions); 00414 } 00415 00416 if (Dim == 3) 00417 { 00418 this->dpsideta_map.resize (n_mapping_shape_functions); 00419 this->d2psidxideta_map.resize (n_mapping_shape_functions); 00420 this->d2psideta2_map.resize (n_mapping_shape_functions); 00421 } 00422 00423 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00424 { 00425 // Allocate space to store the values of the shape functions 00426 // and their first and second derivatives at the quadrature points. 00427 this->psi_map[i].resize (n_qp); 00428 if (Dim > 1) 00429 { 00430 this->dpsidxi_map[i].resize (n_qp); 00431 this->d2psidxi2_map[i].resize (n_qp); 00432 } 00433 if (Dim == 3) 00434 { 00435 this->dpsideta_map[i].resize (n_qp); 00436 this->d2psidxideta_map[i].resize (n_qp); 00437 this->d2psideta2_map[i].resize (n_qp); 00438 } 00439 00440 // Compute the value of shape function i, and its first and 00441 // second derivatives at quadrature point p 00442 // (Lagrange shape functions are used for the mapping) 00443 for (unsigned int p=0; p<n_qp; p++) 00444 { 00445 this->psi_map[i][p] = FE<Dim-1,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00446 if (Dim > 1) 00447 { 00448 this->dpsidxi_map[i][p] = FE<Dim-1,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00449 this->d2psidxi2_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 0, qp[p]); 00450 } 00451 // libMesh::out << "this->d2psidxi2_map["<<i<<"][p]=" << d2psidxi2_map[i][p] << std::endl; 00452 00453 // If we are in 3D, then our sides are 2D faces. 00454 // For the second derivatives, we must also compute the cross 00455 // derivative d^2() / dxi deta 00456 if (Dim == 3) 00457 { 00458 this->dpsideta_map[i][p] = FE<Dim-1,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]); 00459 this->d2psidxideta_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 1, qp[p]); 00460 this->d2psideta2_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 2, qp[p]); 00461 } 00462 } 00463 } 00464 00465 00469 STOP_LOG("init_face_shape_functions()", "FEMap"); 00470 } 00471 00472 template<unsigned int Dim> 00473 void FEMap::init_edge_shape_functions(const std::vector<Point>& qp, 00474 const Elem* edge) 00475 { 00476 libmesh_assert(edge); 00477 00481 START_LOG("init_edge_shape_functions()", "FEMap"); 00482 00483 // The element type and order to use in 00484 // the map 00485 const Order mapping_order (edge->default_order()); 00486 const ElemType mapping_elem_type (edge->type()); 00487 00488 // The number of quadrature points. 00489 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qp.size()); 00490 00491 const unsigned int n_mapping_shape_functions = 00492 FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type, 00493 mapping_order); 00494 00495 // resize the vectors to hold current data 00496 // Psi are the shape functions used for the FE mapping 00497 this->psi_map.resize (n_mapping_shape_functions); 00498 this->dpsidxi_map.resize (n_mapping_shape_functions); 00499 this->d2psidxi2_map.resize (n_mapping_shape_functions); 00500 00501 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00502 { 00503 // Allocate space to store the values of the shape functions 00504 // and their first and second derivatives at the quadrature points. 00505 this->psi_map[i].resize (n_qp); 00506 this->dpsidxi_map[i].resize (n_qp); 00507 this->d2psidxi2_map[i].resize (n_qp); 00508 00509 // Compute the value of shape function i, and its first and 00510 // second derivatives at quadrature point p 00511 // (Lagrange shape functions are used for the mapping) 00512 for (unsigned int p=0; p<n_qp; p++) 00513 { 00514 this->psi_map[i][p] = FE<1,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00515 this->dpsidxi_map[i][p] = FE<1,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00516 this->d2psidxi2_map[i][p] = FE<1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 0, qp[p]); 00517 } 00518 } 00519 00523 STOP_LOG("init_edge_shape_functions()", "FEMap"); 00524 } 00525 00526 00527 00528 void FEMap::compute_face_map(int dim, const std::vector<Real>& qw, 00529 const Elem* side) 00530 { 00531 libmesh_assert(side); 00532 00533 START_LOG("compute_face_map()", "FEMap"); 00534 00535 // The number of quadrature points. 00536 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qw.size()); 00537 00538 switch (dim) 00539 { 00540 case 1: 00541 { 00542 // A 1D finite element, currently assumed to be in 1D space 00543 // This means the boundary is a "0D finite element", a 00544 // NODEELEM. 00545 00546 // Resize the vectors to hold data at the quadrature points 00547 { 00548 this->xyz.resize(n_qp); 00549 normals.resize(n_qp); 00550 00551 this->JxW.resize(n_qp); 00552 } 00553 00554 // If we have no quadrature points, there's nothing else to do 00555 if (!n_qp) 00556 break; 00557 00558 // We need to look back at the full edge to figure out the normal 00559 // vector 00560 const Elem *elem = side->parent(); 00561 libmesh_assert (elem); 00562 if (side->node(0) == elem->node(0)) 00563 normals[0] = Point(-1.); 00564 else 00565 { 00566 libmesh_assert_equal_to (side->node(0), elem->node(1)); 00567 normals[0] = Point(1.); 00568 } 00569 00570 // Calculate x at the point 00571 libmesh_assert_equal_to (this->psi_map.size(), 1); 00572 // In the unlikely event we have multiple quadrature 00573 // points, they'll be in the same place 00574 for (unsigned int p=0; p<n_qp; p++) 00575 { 00576 this->xyz[p].zero(); 00577 this->xyz[p].add_scaled (side->point(0), this->psi_map[0][p]); 00578 normals[p] = normals[0]; 00579 this->JxW[p] = 1.0*qw[p]; 00580 } 00581 00582 // done computing the map 00583 break; 00584 } 00585 00586 case 2: 00587 { 00588 // A 2D finite element living in either 2D or 3D space. 00589 // This means the boundary is a 1D finite element, i.e. 00590 // and EDGE2 or EDGE3. 00591 // Resize the vectors to hold data at the quadrature points 00592 { 00593 this->xyz.resize(n_qp); 00594 this->dxyzdxi_map.resize(n_qp); 00595 this->d2xyzdxi2_map.resize(n_qp); 00596 this->tangents.resize(n_qp); 00597 this->normals.resize(n_qp); 00598 this->curvatures.resize(n_qp); 00599 00600 this->JxW.resize(n_qp); 00601 } 00602 00603 // Clear the entities that will be summed 00604 // Compute the tangent & normal at the quadrature point 00605 for (unsigned int p=0; p<n_qp; p++) 00606 { 00607 this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D 00608 this->xyz[p].zero(); 00609 this->dxyzdxi_map[p].zero(); 00610 this->d2xyzdxi2_map[p].zero(); 00611 } 00612 00613 // compute x, dxdxi at the quadrature points 00614 for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes 00615 { 00616 const Point& side_point = side->point(i); 00617 00618 for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... 00619 { 00620 this->xyz[p].add_scaled (side_point, this->psi_map[i][p]); 00621 this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]); 00622 this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]); 00623 } 00624 } 00625 00626 // Compute the tangent & normal at the quadrature point 00627 for (unsigned int p=0; p<n_qp; p++) 00628 { 00629 // The first tangent comes from just the edge's Jacobian 00630 this->tangents[p][0] = this->dxyzdxi_map[p].unit(); 00631 00632 #if LIBMESH_DIM == 2 00633 // For a 2D element living in 2D, the normal is given directly 00634 // from the entries in the edge Jacobian. 00635 this->normals[p] = (Point(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.)).unit(); 00636 00637 #elif LIBMESH_DIM == 3 00638 // For a 2D element living in 3D, there is a second tangent. 00639 // For the second tangent, we need to refer to the full 00640 // element's (not just the edge's) Jacobian. 00641 const Elem *elem = side->parent(); 00642 libmesh_assert(elem); 00643 00644 // Inverse map xyz[p] to a reference point on the parent... 00645 Point reference_point = FE<2,LAGRANGE>::inverse_map(elem, this->xyz[p]); 00646 00647 // Get dxyz/dxi and dxyz/deta from the parent map. 00648 Point dx_dxi = FE<2,LAGRANGE>::map_xi (elem, reference_point); 00649 Point dx_deta = FE<2,LAGRANGE>::map_eta(elem, reference_point); 00650 00651 // The second tangent vector is formed by crossing these vectors. 00652 tangents[p][1] = dx_dxi.cross(dx_deta).unit(); 00653 00654 // Finally, the normal in this case is given by crossing these 00655 // two tangents. 00656 normals[p] = tangents[p][0].cross(tangents[p][1]).unit(); 00657 #endif 00658 00659 00660 // The curvature is computed via the familiar Frenet formula: 00661 // curvature = [d^2(x) / d (xi)^2] dot [normal] 00662 // For a reference, see: 00663 // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310 00664 // 00665 // Note: The sign convention here is different from the 00666 // 3D case. Concave-upward curves (smiles) have a positive 00667 // curvature. Concave-downward curves (frowns) have a 00668 // negative curvature. Be sure to take that into account! 00669 const Real numerator = this->d2xyzdxi2_map[p] * this->normals[p]; 00670 const Real denominator = this->dxyzdxi_map[p].size_sq(); 00671 libmesh_assert_not_equal_to (denominator, 0); 00672 curvatures[p] = numerator / denominator; 00673 } 00674 00675 // compute the jacobian at the quadrature points 00676 for (unsigned int p=0; p<n_qp; p++) 00677 { 00678 const Real the_jac = this->dxyzdxi_map[p].size(); 00679 00680 libmesh_assert_greater (the_jac, 0.); 00681 00682 this->JxW[p] = the_jac*qw[p]; 00683 } 00684 00685 // done computing the map 00686 break; 00687 } 00688 00689 00690 00691 case 3: 00692 { 00693 // A 3D finite element living in 3D space. 00694 // Resize the vectors to hold data at the quadrature points 00695 { 00696 this->xyz.resize(n_qp); 00697 this->dxyzdxi_map.resize(n_qp); 00698 this->dxyzdeta_map.resize(n_qp); 00699 this->d2xyzdxi2_map.resize(n_qp); 00700 this->d2xyzdxideta_map.resize(n_qp); 00701 this->d2xyzdeta2_map.resize(n_qp); 00702 this->tangents.resize(n_qp); 00703 this->normals.resize(n_qp); 00704 this->curvatures.resize(n_qp); 00705 00706 this->JxW.resize(n_qp); 00707 } 00708 00709 // Clear the entities that will be summed 00710 for (unsigned int p=0; p<n_qp; p++) 00711 { 00712 this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D 00713 this->xyz[p].zero(); 00714 this->dxyzdxi_map[p].zero(); 00715 this->dxyzdeta_map[p].zero(); 00716 this->d2xyzdxi2_map[p].zero(); 00717 this->d2xyzdxideta_map[p].zero(); 00718 this->d2xyzdeta2_map[p].zero(); 00719 } 00720 00721 // compute x, dxdxi at the quadrature points 00722 for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes 00723 { 00724 const Point& side_point = side->point(i); 00725 00726 for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... 00727 { 00728 this->xyz[p].add_scaled (side_point, this->psi_map[i][p]); 00729 this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]); 00730 this->dxyzdeta_map[p].add_scaled(side_point, this->dpsideta_map[i][p]); 00731 this->d2xyzdxi2_map[p].add_scaled (side_point, this->d2psidxi2_map[i][p]); 00732 this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]); 00733 this->d2xyzdeta2_map[p].add_scaled (side_point, this->d2psideta2_map[i][p]); 00734 } 00735 } 00736 00737 // Compute the tangents, normal, and curvature at the quadrature point 00738 for (unsigned int p=0; p<n_qp; p++) 00739 { 00740 const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]); 00741 this->normals[p] = n.unit(); 00742 this->tangents[p][0] = this->dxyzdxi_map[p].unit(); 00743 this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit(); 00744 00745 // Compute curvature using the typical nomenclature 00746 // of the first and second fundamental forms. 00747 // For reference, see: 00748 // 1) http://mathworld.wolfram.com/MeanCurvature.html 00749 // (note -- they are using inward normal) 00750 // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill 00751 const Real L = -this->d2xyzdxi2_map[p] * this->normals[p]; 00752 const Real M = -this->d2xyzdxideta_map[p] * this->normals[p]; 00753 const Real N = -this->d2xyzdeta2_map[p] * this->normals[p]; 00754 const Real E = this->dxyzdxi_map[p].size_sq(); 00755 const Real F = this->dxyzdxi_map[p] * this->dxyzdeta_map[p]; 00756 const Real G = this->dxyzdeta_map[p].size_sq(); 00757 00758 const Real numerator = E*N -2.*F*M + G*L; 00759 const Real denominator = E*G - F*F; 00760 libmesh_assert_not_equal_to (denominator, 0.); 00761 curvatures[p] = 0.5*numerator/denominator; 00762 } 00763 00764 // compute the jacobian at the quadrature points, see 00765 // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm 00766 for (unsigned int p=0; p<n_qp; p++) 00767 { 00768 const Real g11 = (dxdxi_map(p)*dxdxi_map(p) + 00769 dydxi_map(p)*dydxi_map(p) + 00770 dzdxi_map(p)*dzdxi_map(p)); 00771 00772 const Real g12 = (dxdxi_map(p)*dxdeta_map(p) + 00773 dydxi_map(p)*dydeta_map(p) + 00774 dzdxi_map(p)*dzdeta_map(p)); 00775 00776 const Real g21 = g12; 00777 00778 const Real g22 = (dxdeta_map(p)*dxdeta_map(p) + 00779 dydeta_map(p)*dydeta_map(p) + 00780 dzdeta_map(p)*dzdeta_map(p)); 00781 00782 00783 const Real the_jac = std::sqrt(g11*g22 - g12*g21); 00784 00785 libmesh_assert_greater (the_jac, 0.); 00786 00787 this->JxW[p] = the_jac*qw[p]; 00788 } 00789 00790 // done computing the map 00791 break; 00792 } 00793 00794 00795 default: 00796 libmesh_error(); 00797 00798 } 00799 STOP_LOG("compute_face_map()", "FEMap"); 00800 } 00801 00802 00803 00804 00805 void FEMap::compute_edge_map(int dim, 00806 const std::vector<Real>& qw, 00807 const Elem* edge) 00808 { 00809 libmesh_assert(edge); 00810 00811 if (dim == 2) 00812 { 00813 // A 2D finite element living in either 2D or 3D space. 00814 // The edges here are the sides of the element, so the 00815 // (misnamed) compute_face_map function does what we want 00816 this->compute_face_map(dim, qw, edge); 00817 return; 00818 } 00819 00820 libmesh_assert_equal_to (dim, 3); // 1D is unnecessary and currently unsupported 00821 00822 START_LOG("compute_edge_map()", "FEMap"); 00823 00824 // The number of quadrature points. 00825 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qw.size()); 00826 00827 // Resize the vectors to hold data at the quadrature points 00828 this->xyz.resize(n_qp); 00829 this->dxyzdxi_map.resize(n_qp); 00830 this->dxyzdeta_map.resize(n_qp); 00831 this->d2xyzdxi2_map.resize(n_qp); 00832 this->d2xyzdxideta_map.resize(n_qp); 00833 this->d2xyzdeta2_map.resize(n_qp); 00834 this->tangents.resize(n_qp); 00835 this->normals.resize(n_qp); 00836 this->curvatures.resize(n_qp); 00837 00838 this->JxW.resize(n_qp); 00839 00840 // Clear the entities that will be summed 00841 for (unsigned int p=0; p<n_qp; p++) 00842 { 00843 this->tangents[p].resize(1); 00844 this->xyz[p].zero(); 00845 this->dxyzdxi_map[p].zero(); 00846 this->dxyzdeta_map[p].zero(); 00847 this->d2xyzdxi2_map[p].zero(); 00848 this->d2xyzdxideta_map[p].zero(); 00849 this->d2xyzdeta2_map[p].zero(); 00850 } 00851 00852 // compute x, dxdxi at the quadrature points 00853 for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes 00854 { 00855 const Point& edge_point = edge->point(i); 00856 00857 for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... 00858 { 00859 this->xyz[p].add_scaled (edge_point, this->psi_map[i][p]); 00860 this->dxyzdxi_map[p].add_scaled (edge_point, this->dpsidxi_map[i][p]); 00861 this->d2xyzdxi2_map[p].add_scaled (edge_point, this->d2psidxi2_map[i][p]); 00862 } 00863 } 00864 00865 // Compute the tangents at the quadrature point 00866 // FIXME: normals (plural!) and curvatures are uncalculated 00867 for (unsigned int p=0; p<n_qp; p++) 00868 { 00869 const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]); 00870 this->tangents[p][0] = this->dxyzdxi_map[p].unit(); 00871 00872 // compute the jacobian at the quadrature points 00873 const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) + 00874 this->dydxi_map(p)*this->dydxi_map(p) + 00875 this->dzdxi_map(p)*this->dzdxi_map(p)); 00876 00877 libmesh_assert_greater (the_jac, 0.); 00878 00879 this->JxW[p] = the_jac*qw[p]; 00880 } 00881 00882 STOP_LOG("compute_edge_map()", "FEMap"); 00883 } 00884 00885 00886 // Explicit FEMap Instantiations 00887 FACE_EDGE_SHAPE_ERROR(0,init_face_shape_functions) 00888 template void FEMap::init_face_shape_functions<1>(const std::vector<Point>&,const Elem*); 00889 template void FEMap::init_face_shape_functions<2>(const std::vector<Point>&,const Elem*); 00890 template void FEMap::init_face_shape_functions<3>(const std::vector<Point>&,const Elem*); 00891 00892 FACE_EDGE_SHAPE_ERROR(0,init_edge_shape_functions) 00893 template void FEMap::init_edge_shape_functions<1>(const std::vector<Point>&, const Elem*); 00894 template void FEMap::init_edge_shape_functions<2>(const std::vector<Point>&, const Elem*); 00895 template void FEMap::init_edge_shape_functions<3>(const std::vector<Point>&, const Elem*); 00896 00897 //-------------------------------------------------------------- 00898 // Explicit FE instantiations 00899 template void FE<1,LAGRANGE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00900 template void FE<1,LAGRANGE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00901 template void FE<1,LAGRANGE_VEC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00902 template void FE<1,LAGRANGE_VEC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00903 template void FE<1,L2_LAGRANGE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00904 template void FE<1,L2_LAGRANGE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00905 template void FE<1,HIERARCHIC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00906 template void FE<1,HIERARCHIC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00907 template void FE<1,L2_HIERARCHIC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00908 template void FE<1,L2_HIERARCHIC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00909 template void FE<1,CLOUGH>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00910 template void FE<1,CLOUGH>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00911 template void FE<1,HERMITE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00912 template void FE<1,HERMITE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00913 template void FE<1,MONOMIAL>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00914 template void FE<1,MONOMIAL>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00915 template void FE<1,SCALAR>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00916 template void FE<1,SCALAR>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00917 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00918 template void FE<1,BERNSTEIN>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00919 template void FE<1,BERNSTEIN>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00920 template void FE<1,SZABAB>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00921 template void FE<1,SZABAB>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00922 #endif 00923 template void FE<1,XYZ>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00924 template void FE<1,XYZ>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00925 00926 template void FE<2,LAGRANGE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00927 template void FE<2,LAGRANGE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00928 template void FE<2,LAGRANGE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00929 template void FE<2,LAGRANGE_VEC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00930 template void FE<2,LAGRANGE_VEC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00931 template void FE<2,LAGRANGE_VEC>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00932 template void FE<2,L2_LAGRANGE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00933 template void FE<2,L2_LAGRANGE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00934 template void FE<2,L2_LAGRANGE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00935 template void FE<2,HIERARCHIC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00936 template void FE<2,HIERARCHIC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00937 template void FE<2,HIERARCHIC>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00938 template void FE<2,L2_HIERARCHIC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00939 template void FE<2,L2_HIERARCHIC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00940 template void FE<2,L2_HIERARCHIC>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00941 template void FE<2,CLOUGH>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00942 template void FE<2,CLOUGH>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00943 template void FE<2,CLOUGH>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00944 template void FE<2,HERMITE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00945 template void FE<2,HERMITE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00946 template void FE<2,HERMITE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00947 template void FE<2,MONOMIAL>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00948 template void FE<2,MONOMIAL>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00949 template void FE<2,MONOMIAL>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00950 template void FE<2,SCALAR>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00951 template void FE<2,SCALAR>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00952 template void FE<2,SCALAR>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00953 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00954 template void FE<2,BERNSTEIN>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00955 template void FE<2,BERNSTEIN>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00956 template void FE<2,BERNSTEIN>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00957 template void FE<2,SZABAB>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00958 template void FE<2,SZABAB>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00959 template void FE<2,SZABAB>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00960 #endif 00961 template void FE<2,XYZ>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00962 template void FE<2,XYZ>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00963 template void FE<2,XYZ>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00964 template void FE<2,NEDELEC_ONE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00965 template void FE<2,NEDELEC_ONE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00966 template void FE<2,NEDELEC_ONE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00967 00968 // Intel 9.1 complained it needed this in devel mode. 00969 //template void FE<2,XYZ>::init_face_shape_functions(const std::vector<Point>&, const Elem*); 00970 00971 template void FE<3,LAGRANGE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00972 template void FE<3,LAGRANGE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00973 template void FE<3,LAGRANGE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00974 template void FE<3,LAGRANGE_VEC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00975 template void FE<3,LAGRANGE_VEC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00976 template void FE<3,LAGRANGE_VEC>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00977 template void FE<3,L2_LAGRANGE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00978 template void FE<3,L2_LAGRANGE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00979 template void FE<3,L2_LAGRANGE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00980 template void FE<3,HIERARCHIC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00981 template void FE<3,HIERARCHIC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00982 template void FE<3,HIERARCHIC>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00983 template void FE<3,L2_HIERARCHIC>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00984 template void FE<3,L2_HIERARCHIC>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00985 template void FE<3,L2_HIERARCHIC>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00986 template void FE<3,CLOUGH>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00987 template void FE<3,CLOUGH>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00988 template void FE<3,CLOUGH>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00989 template void FE<3,HERMITE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00990 template void FE<3,HERMITE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00991 template void FE<3,HERMITE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00992 template void FE<3,MONOMIAL>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00993 template void FE<3,MONOMIAL>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00994 template void FE<3,MONOMIAL>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00995 template void FE<3,SCALAR>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00996 template void FE<3,SCALAR>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 00997 template void FE<3,SCALAR>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 00998 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00999 template void FE<3,BERNSTEIN>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01000 template void FE<3,BERNSTEIN>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 01001 template void FE<3,BERNSTEIN>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01002 template void FE<3,SZABAB>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01003 template void FE<3,SZABAB>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 01004 template void FE<3,SZABAB>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01005 #endif 01006 template void FE<3,XYZ>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01007 template void FE<3,XYZ>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 01008 template void FE<3,XYZ>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01009 template void FE<3,NEDELEC_ONE>::reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01010 template void FE<3,NEDELEC_ONE>::side_map(Elem const*, Elem const*, const unsigned int, const std::vector<Point>&, std::vector<Point>&); 01011 template void FE<3,NEDELEC_ONE>::edge_reinit(Elem const*, unsigned int, Real, const std::vector<Point>* const, const std::vector<Real>* const); 01012 01013 // Intel 9.1 complained it needed this in devel mode. 01014 //template void FE<3,XYZ>::init_face_shape_functions(const std::vector<Point>&, const Elem*); 01015 01016 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: