$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 // Local includes 00021 #include "libmesh/dof_map.h" 00022 #include "libmesh/fe.h" 00023 #include "libmesh/fe_interface.h" 00024 #include "libmesh/elem.h" 00025 #include "libmesh/threads.h" 00026 #include "libmesh/tensor_value.h" 00027 #include "libmesh/string_to_enum.h" 00028 00029 namespace libMesh 00030 { 00031 00032 // ------------------------------------------------------------ 00033 // Nedelec first kind specific implementations 00034 00035 00036 // Anonymous namespace for local helper functions 00037 namespace { 00038 void nedelec_one_nodal_soln(const Elem* elem, 00039 const Order order, 00040 const std::vector<Number>& elem_soln, 00041 const int dim, 00042 std::vector<Number>& nodal_soln) 00043 { 00044 const unsigned int n_nodes = elem->n_nodes(); 00045 const ElemType elem_type = elem->type(); 00046 00047 const Order totalorder = static_cast<Order>(order+elem->p_level()); 00048 00049 nodal_soln.resize(n_nodes*dim); 00050 00051 FEType fe_type(totalorder, NEDELEC_ONE); 00052 00053 switch (totalorder) 00054 { 00055 case FIRST: 00056 { 00057 switch (elem_type) 00058 { 00059 case TRI6: 00060 { 00061 libmesh_assert_equal_to (elem_soln.size(), 3); 00062 libmesh_assert_equal_to (nodal_soln.size(), 6*2); 00063 break; 00064 } 00065 case QUAD8: 00066 case QUAD9: 00067 { 00068 libmesh_assert_equal_to (elem_soln.size(), 4); 00069 00070 if (elem_type == QUAD8) 00071 libmesh_assert_equal_to (nodal_soln.size(), 8*2); 00072 else 00073 libmesh_assert_equal_to (nodal_soln.size(), 9*2); 00074 break; 00075 } 00076 case TET10: 00077 { 00078 libmesh_assert_equal_to (elem_soln.size(), 6); 00079 libmesh_assert_equal_to (nodal_soln.size(), 10*3); 00080 00081 libmesh_not_implemented(); 00082 00083 break; 00084 } 00085 00086 00087 case HEX20: 00088 case HEX27: 00089 { 00090 libmesh_assert_equal_to (elem_soln.size(), 12); 00091 00092 if (elem_type == HEX20) 00093 libmesh_assert_equal_to (nodal_soln.size(), 20*3); 00094 else 00095 libmesh_assert_equal_to (nodal_soln.size(), 27*3); 00096 00097 break; 00098 } 00099 00100 default: 00101 libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(elem_type) << " selected for NEDELEC_ONE FE family!"); 00102 00103 } // switch(elem_type) 00104 00105 const unsigned int n_sf = 00106 FEInterface::n_shape_functions(dim, fe_type, elem_type); 00107 00108 std::vector<Point> refspace_nodes; 00109 FEVectorBase::get_refspace_nodes(elem_type,refspace_nodes); 00110 libmesh_assert_equal_to (refspace_nodes.size(), n_nodes); 00111 00112 00113 // Need to create new fe object so the shape function as the FETransformation 00114 // applied to it. 00115 UniquePtr<FEVectorBase> vis_fe = FEVectorBase::build(dim,fe_type); 00116 00117 const std::vector<std::vector<RealGradient> >& vis_phi = vis_fe->get_phi(); 00118 00119 vis_fe->reinit(elem,&refspace_nodes); 00120 00121 for( unsigned int n = 0; n < n_nodes; n++ ) 00122 { 00123 libmesh_assert_equal_to (elem_soln.size(), n_sf); 00124 00125 // Zero before summation 00126 for( int d = 0; d < dim; d++ ) 00127 { 00128 nodal_soln[dim*n+d] = 0; 00129 } 00130 00131 // u = Sum (u_i phi_i) 00132 for (unsigned int i=0; i<n_sf; i++) 00133 { 00134 for( int d = 0; d < dim; d++ ) 00135 { 00136 nodal_soln[dim*n+d] += elem_soln[i]*(vis_phi[i][n](d)); 00137 } 00138 } 00139 } 00140 00141 return; 00142 } // case FIRST 00143 00144 default: 00145 libmesh_error_msg("ERROR: Invalid total order " << Utility::enum_to_string(totalorder) << " selected for NEDELEC_ONE FE family!"); 00146 00147 }//switch (totalorder) 00148 00149 return; 00150 } // nedelec_one_nodal_soln 00151 00152 00153 unsigned int nedelec_one_n_dofs(const ElemType t, const Order o) 00154 { 00155 switch (o) 00156 { 00157 case FIRST: 00158 { 00159 switch (t) 00160 { 00161 case TRI6: 00162 return 3; 00163 00164 case QUAD8: 00165 case QUAD9: 00166 return 4; 00167 00168 case TET10: 00169 return 6; 00170 00171 case HEX20: 00172 case HEX27: 00173 return 12; 00174 00175 case INVALID_ELEM: 00176 return 0; 00177 00178 default: 00179 libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!"); 00180 } 00181 } 00182 00183 default: 00184 libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for NEDELEC_ONE FE family!"); 00185 } 00186 00187 libmesh_error_msg("We'll never get here!"); 00188 return 0; 00189 } 00190 00191 00192 00193 00194 unsigned int nedelec_one_n_dofs_at_node(const ElemType t, 00195 const Order o, 00196 const unsigned int n) 00197 { 00198 switch (o) 00199 { 00200 case FIRST: 00201 { 00202 switch (t) 00203 { 00204 case TRI6: 00205 { 00206 switch (n) 00207 { 00208 case 0: 00209 case 1: 00210 case 2: 00211 return 0; 00212 case 3: 00213 case 4: 00214 case 5: 00215 return 1; 00216 00217 default: 00218 libmesh_error_msg("ERROR: Invalid node ID " << n); 00219 } 00220 } 00221 case QUAD8: 00222 { 00223 switch (n) 00224 { 00225 case 0: 00226 case 1: 00227 case 2: 00228 case 3: 00229 return 0; 00230 case 4: 00231 case 5: 00232 case 6: 00233 case 7: 00234 return 1; 00235 00236 default: 00237 libmesh_error_msg("ERROR: Invalid node ID " << n); 00238 } 00239 } 00240 case QUAD9: 00241 { 00242 switch (n) 00243 { 00244 case 0: 00245 case 1: 00246 case 2: 00247 case 3: 00248 case 8: 00249 return 0; 00250 case 4: 00251 case 5: 00252 case 6: 00253 case 7: 00254 return 1; 00255 00256 default: 00257 libmesh_error_msg("ERROR: Invalid node ID " << n); 00258 } 00259 } 00260 case TET10: 00261 { 00262 switch (n) 00263 { 00264 case 0: 00265 case 1: 00266 case 2: 00267 case 3: 00268 return 0; 00269 case 4: 00270 case 5: 00271 case 6: 00272 case 7: 00273 case 8: 00274 case 9: 00275 return 1; 00276 00277 default: 00278 libmesh_error_msg("ERROR: Invalid node ID " << n); 00279 } 00280 } 00281 00282 case HEX20: 00283 { 00284 switch (n) 00285 { 00286 case 0: 00287 case 1: 00288 case 2: 00289 case 3: 00290 case 4: 00291 case 5: 00292 case 6: 00293 case 7: 00294 return 0; 00295 case 8: 00296 case 9: 00297 case 10: 00298 case 11: 00299 case 12: 00300 case 13: 00301 case 14: 00302 case 15: 00303 case 16: 00304 case 17: 00305 case 18: 00306 case 19: 00307 return 1; 00308 00309 default: 00310 libmesh_error_msg("ERROR: Invalid node ID " << n); 00311 } 00312 } 00313 case HEX27: 00314 { 00315 switch (n) 00316 { 00317 case 0: 00318 case 1: 00319 case 2: 00320 case 3: 00321 case 4: 00322 case 5: 00323 case 6: 00324 case 7: 00325 case 20: 00326 case 21: 00327 case 22: 00328 case 23: 00329 case 24: 00330 case 25: 00331 case 26: 00332 return 0; 00333 case 8: 00334 case 9: 00335 case 10: 00336 case 11: 00337 case 12: 00338 case 13: 00339 case 14: 00340 case 15: 00341 case 16: 00342 case 17: 00343 case 18: 00344 case 19: 00345 return 1; 00346 00347 default: 00348 libmesh_error_msg("ERROR: Invalid node ID " << n); 00349 } 00350 } 00351 00352 case INVALID_ELEM: 00353 return 0; 00354 00355 default: 00356 libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!"); 00357 } 00358 } 00359 00360 default: 00361 libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for NEDELEC_ONE FE family!"); 00362 } 00363 00364 libmesh_error_msg("We'll never get here!"); 00365 return 0; 00366 } 00367 00368 00369 00370 #ifdef LIBMESH_ENABLE_AMR 00371 void nedelec_one_compute_constraints (DofConstraints &/*constraints*/, 00372 DofMap &/*dof_map*/, 00373 const unsigned int /*variable_number*/, 00374 const Elem* libmesh_dbg_var(elem), 00375 const unsigned Dim) 00376 { 00377 // Only constrain elements in 2,3D. 00378 if (Dim == 1) 00379 return; 00380 00381 libmesh_assert(elem); 00382 00383 libmesh_not_implemented(); 00384 00385 /* 00386 // Only constrain active and ancestor elements 00387 if (elem->subactive()) 00388 return; 00389 00390 FEType fe_type = dof_map.variable_type(variable_number); 00391 fe_type.order = static_cast<Order>(fe_type.order + elem->p_level()); 00392 00393 std::vector<unsigned int> my_dof_indices, parent_dof_indices; 00394 00395 // Look at the element faces. Check to see if we need to 00396 // build constraints. 00397 for (unsigned int s=0; s<elem->n_sides(); s++) 00398 if (elem->neighbor(s) != NULL) 00399 if (elem->neighbor(s)->level() < elem->level()) // constrain dofs shared between 00400 { // this element and ones coarser 00401 // than this element. 00402 // Get pointers to the elements of interest and its parent. 00403 const Elem* parent = elem->parent(); 00404 00405 // This can't happen... Only level-0 elements have NULL 00406 // parents, and no level-0 elements can be at a higher 00407 // level than their neighbors! 00408 libmesh_assert(parent); 00409 00410 const UniquePtr<Elem> my_side (elem->build_side(s)); 00411 const UniquePtr<Elem> parent_side (parent->build_side(s)); 00412 00413 // This function gets called element-by-element, so there 00414 // will be a lot of memory allocation going on. We can 00415 // at least minimize this for the case of the dof indices 00416 // by efficiently preallocating the requisite storage. 00417 my_dof_indices.reserve (my_side->n_nodes()); 00418 parent_dof_indices.reserve (parent_side->n_nodes()); 00419 00420 dof_map.dof_indices (my_side.get(), my_dof_indices, 00421 variable_number); 00422 dof_map.dof_indices (parent_side.get(), parent_dof_indices, 00423 variable_number); 00424 00425 for (unsigned int my_dof=0; 00426 my_dof<FEInterface::n_dofs(Dim-1, fe_type, my_side->type()); 00427 my_dof++) 00428 { 00429 libmesh_assert_less (my_dof, my_side->n_nodes()); 00430 00431 // My global dof index. 00432 const unsigned int my_dof_g = my_dof_indices[my_dof]; 00433 00434 // The support point of the DOF 00435 const Point& support_point = my_side->point(my_dof); 00436 00437 // Figure out where my node lies on their reference element. 00438 const Point mapped_point = FEInterface::inverse_map(Dim-1, fe_type, 00439 parent_side.get(), 00440 support_point); 00441 00442 // Compute the parent's side shape function values. 00443 for (unsigned int their_dof=0; 00444 their_dof<FEInterface::n_dofs(Dim-1, fe_type, parent_side->type()); 00445 their_dof++) 00446 { 00447 libmesh_assert_less (their_dof, parent_side->n_nodes()); 00448 00449 // Their global dof index. 00450 const unsigned int their_dof_g = 00451 parent_dof_indices[their_dof]; 00452 00453 const Real their_dof_value = FEInterface::shape(Dim-1, 00454 fe_type, 00455 parent_side->type(), 00456 their_dof, 00457 mapped_point); 00458 00459 // Only add non-zero and non-identity values 00460 // for Lagrange basis functions. 00461 if ((std::abs(their_dof_value) > 1.e-5) && 00462 (std::abs(their_dof_value) < .999)) 00463 { 00464 // since we may be running this method concurrently 00465 // on multiple threads we need to acquire a lock 00466 // before modifying the shared constraint_row object. 00467 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00468 00469 // A reference to the constraint row. 00470 DofConstraintRow& constraint_row = constraints[my_dof_g].first; 00471 00472 constraint_row.insert(std::make_pair (their_dof_g, 00473 their_dof_value)); 00474 } 00475 #ifdef DEBUG 00476 // Protect for the case u_i = 0.999 u_j, 00477 // in which case i better equal j. 00478 else if (their_dof_value >= .999) 00479 libmesh_assert_equal_to (my_dof_g, their_dof_g); 00480 #endif 00481 } 00482 } 00483 } 00484 */ 00485 } // nedelec_one_compute_constrants() 00486 #endif // #ifdef LIBMESH_ENABLE_AMR 00487 00488 } // anonymous namespace 00489 00490 #define NEDELEC_LOW_D_ERROR_MESSAGE \ 00491 libmesh_error_msg("ERROR: This method makes no sense for low-D elements!"); 00492 00493 00494 // Do full-specialization for every dimension, instead 00495 // of explicit instantiation at the end of this file. 00496 template <> 00497 void FE<0,NEDELEC_ONE>::nodal_soln(const Elem*, 00498 const Order, 00499 const std::vector<Number>&, 00500 std::vector<Number>&) 00501 { NEDELEC_LOW_D_ERROR_MESSAGE } 00502 00503 template <> 00504 void FE<1,NEDELEC_ONE>::nodal_soln(const Elem*, 00505 const Order, 00506 const std::vector<Number>&, 00507 std::vector<Number>&) 00508 { NEDELEC_LOW_D_ERROR_MESSAGE } 00509 00510 template <> 00511 void FE<2,NEDELEC_ONE>::nodal_soln(const Elem* elem, 00512 const Order order, 00513 const std::vector<Number>& elem_soln, 00514 std::vector<Number>& nodal_soln) 00515 { nedelec_one_nodal_soln(elem, order, elem_soln, 2 /*dim*/, nodal_soln); } 00516 00517 template <> 00518 void FE<3,NEDELEC_ONE>::nodal_soln(const Elem* elem, 00519 const Order order, 00520 const std::vector<Number>& elem_soln, 00521 std::vector<Number>& nodal_soln) 00522 { nedelec_one_nodal_soln(elem, order, elem_soln, 3 /*dim*/, nodal_soln); } 00523 00524 00525 // Do full-specialization for every dimension, instead 00526 // of explicit instantiation at the end of this function. 00527 // This could be macro-ified. 00528 template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE } 00529 template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE } 00530 template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs(const ElemType t, const Order o) { return nedelec_one_n_dofs(t, o); } 00531 template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs(const ElemType t, const Order o) { return nedelec_one_n_dofs(t, o); } 00532 00533 00534 // Do full-specialization for every dimension, instead 00535 // of explicit instantiation at the end of this function. 00536 template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { NEDELEC_LOW_D_ERROR_MESSAGE } 00537 template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { NEDELEC_LOW_D_ERROR_MESSAGE } 00538 template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return nedelec_one_n_dofs_at_node(t, o, n); } 00539 template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return nedelec_one_n_dofs_at_node(t, o, n); } 00540 00541 00542 // Nedelec first type elements have no dofs per element 00543 // FIXME: Only for first order! 00544 template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs_per_elem(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE } 00545 template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs_per_elem(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE } 00546 template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs_per_elem(const ElemType, const Order) { return 0; } 00547 template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs_per_elem(const ElemType, const Order) { return 0; } 00548 00549 // Nedelec first type FEMs are always tangentially continuous 00550 template <> FEContinuity FE<0,NEDELEC_ONE>::get_continuity() const { NEDELEC_LOW_D_ERROR_MESSAGE } 00551 template <> FEContinuity FE<1,NEDELEC_ONE>::get_continuity() const { NEDELEC_LOW_D_ERROR_MESSAGE } 00552 template <> FEContinuity FE<2,NEDELEC_ONE>::get_continuity() const { return H_CURL; } 00553 template <> FEContinuity FE<3,NEDELEC_ONE>::get_continuity() const { return H_CURL; } 00554 00555 // Nedelec first type FEMs are not hierarchic 00556 template <> bool FE<0,NEDELEC_ONE>::is_hierarchic() const { NEDELEC_LOW_D_ERROR_MESSAGE } 00557 template <> bool FE<1,NEDELEC_ONE>::is_hierarchic() const { NEDELEC_LOW_D_ERROR_MESSAGE } 00558 template <> bool FE<2,NEDELEC_ONE>::is_hierarchic() const { return false; } 00559 template <> bool FE<3,NEDELEC_ONE>::is_hierarchic() const { return false; } 00560 00561 // Nedelec first type FEM shapes always need to be reinit'ed (because of orientation dependence) 00562 template <> bool FE<0,NEDELEC_ONE>::shapes_need_reinit() const { NEDELEC_LOW_D_ERROR_MESSAGE } 00563 template <> bool FE<1,NEDELEC_ONE>::shapes_need_reinit() const { NEDELEC_LOW_D_ERROR_MESSAGE } 00564 template <> bool FE<2,NEDELEC_ONE>::shapes_need_reinit() const { return true; } 00565 template <> bool FE<3,NEDELEC_ONE>::shapes_need_reinit() const { return true; } 00566 00567 #ifdef LIBMESH_ENABLE_AMR 00568 template <> 00569 void FE<0,NEDELEC_ONE>::compute_constraints (DofConstraints &, 00570 DofMap &, 00571 const unsigned int, 00572 const Elem*) 00573 { NEDELEC_LOW_D_ERROR_MESSAGE } 00574 00575 template <> 00576 void FE<1,NEDELEC_ONE>::compute_constraints (DofConstraints &, 00577 DofMap &, 00578 const unsigned int, 00579 const Elem*) 00580 { NEDELEC_LOW_D_ERROR_MESSAGE } 00581 00582 template <> 00583 void FE<2,NEDELEC_ONE>::compute_constraints (DofConstraints &constraints, 00584 DofMap &dof_map, 00585 const unsigned int variable_number, 00586 const Elem* elem) 00587 { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/2); } 00588 00589 template <> 00590 void FE<3,NEDELEC_ONE>::compute_constraints (DofConstraints &constraints, 00591 DofMap &dof_map, 00592 const unsigned int variable_number, 00593 const Elem* elem) 00594 { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/3); } 00595 #endif // LIBMESH_ENABLE_AMR 00596 00597 // Specialize useless shape function methods 00598 template <> 00599 RealGradient FE<0,NEDELEC_ONE>::shape(const ElemType, const Order,const unsigned int,const Point&) 00600 { NEDELEC_LOW_D_ERROR_MESSAGE } 00601 template <> 00602 RealGradient FE<0,NEDELEC_ONE>::shape(const Elem*,const Order,const unsigned int,const Point&) 00603 { NEDELEC_LOW_D_ERROR_MESSAGE } 00604 template <> 00605 RealGradient FE<0,NEDELEC_ONE>::shape_deriv(const ElemType, const Order,const unsigned int, 00606 const unsigned int,const Point&) 00607 { NEDELEC_LOW_D_ERROR_MESSAGE } 00608 template <> 00609 RealGradient FE<0,NEDELEC_ONE>::shape_deriv(const Elem*,const Order,const unsigned int, 00610 const unsigned int,const Point&) 00611 { NEDELEC_LOW_D_ERROR_MESSAGE } 00612 template <> 00613 RealGradient FE<0,NEDELEC_ONE>::shape_second_deriv(const ElemType, const Order,const unsigned int, 00614 const unsigned int,const Point&) 00615 { NEDELEC_LOW_D_ERROR_MESSAGE } 00616 template <> 00617 RealGradient FE<0,NEDELEC_ONE>::shape_second_deriv(const Elem*,const Order,const unsigned int, 00618 const unsigned int,const Point&) 00619 { NEDELEC_LOW_D_ERROR_MESSAGE } 00620 00621 template <> 00622 RealGradient FE<1,NEDELEC_ONE>::shape(const ElemType, const Order,const unsigned int,const Point&) 00623 { NEDELEC_LOW_D_ERROR_MESSAGE } 00624 template <> 00625 RealGradient FE<1,NEDELEC_ONE>::shape(const Elem*,const Order,const unsigned int,const Point&) 00626 { NEDELEC_LOW_D_ERROR_MESSAGE } 00627 template <> 00628 RealGradient FE<1,NEDELEC_ONE>::shape_deriv(const ElemType, const Order,const unsigned int, 00629 const unsigned int,const Point&) 00630 { NEDELEC_LOW_D_ERROR_MESSAGE } 00631 template <> 00632 RealGradient FE<1,NEDELEC_ONE>::shape_deriv(const Elem*,const Order,const unsigned int, 00633 const unsigned int,const Point&) 00634 { NEDELEC_LOW_D_ERROR_MESSAGE } 00635 template <> 00636 RealGradient FE<1,NEDELEC_ONE>::shape_second_deriv(const ElemType, const Order,const unsigned int, 00637 const unsigned int,const Point&) 00638 { NEDELEC_LOW_D_ERROR_MESSAGE } 00639 template <> 00640 RealGradient FE<1,NEDELEC_ONE>::shape_second_deriv(const Elem*,const Order,const unsigned int, 00641 const unsigned int,const Point&) 00642 { NEDELEC_LOW_D_ERROR_MESSAGE } 00643 00644 } // namespace libMesh