fe_base.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 // Local includes 00021 #include "libmesh/fe.h" 00022 #include "libmesh/inf_fe.h" 00023 #include "libmesh/libmesh_logging.h" 00024 // For projection code: 00025 #include "libmesh/boundary_info.h" 00026 #include "libmesh/mesh_base.h" 00027 #include "libmesh/dense_matrix.h" 00028 #include "libmesh/dense_vector.h" 00029 #include "libmesh/dof_map.h" 00030 #include "libmesh/elem.h" 00031 #include "libmesh/fe_interface.h" 00032 #include "libmesh/numeric_vector.h" 00033 #include "libmesh/periodic_boundary_base.h" 00034 #include "libmesh/periodic_boundaries.h" 00035 #include "libmesh/quadrature.h" 00036 #include "libmesh/quadrature_gauss.h" 00037 #include "libmesh/tensor_value.h" 00038 #include "libmesh/threads.h" 00039 00040 // Anonymous namespace, for a helper function for periodic boundary 00041 // constraint calculations 00042 namespace 00043 { 00044 using namespace libMesh; 00045 00046 // Find the "primary" element around a boundary point: 00047 const Elem* primary_boundary_point_neighbor 00048 (const Elem* elem, 00049 const Point& p, 00050 const BoundaryInfo& boundary_info, 00051 const std::set<boundary_id_type>& boundary_ids) 00052 { 00053 // If we don't find a better alternative, the user will have 00054 // provided the primary element 00055 const Elem *primary = elem; 00056 00057 std::set<const Elem*> point_neighbors; 00058 elem->find_point_neighbors(p, point_neighbors); 00059 for (std::set<const Elem*>::const_iterator point_neighbors_iter = 00060 point_neighbors.begin(); 00061 point_neighbors_iter != point_neighbors.end(); ++point_neighbors_iter) 00062 { 00063 const Elem* pt_neighbor = *point_neighbors_iter; 00064 00065 // If this point neighbor isn't at least 00066 // as coarse as the current primary elem, or if it is at 00067 // the same level but has a lower id, then 00068 // we won't defer to it. 00069 if ((pt_neighbor->level() > primary->level()) || 00070 (pt_neighbor->level() == primary->level() && 00071 pt_neighbor->id() < primary->id())) 00072 continue; 00073 00074 // Otherwise, we will defer to the point neighbor, but only if 00075 // one of its sides is on a relevant boundary and that side 00076 // contains this vertex 00077 bool vertex_on_periodic_side = false; 00078 for (unsigned int ns = 0; 00079 ns != pt_neighbor->n_sides(); ++ns) 00080 { 00081 const std::vector<boundary_id_type> bc_ids = 00082 boundary_info.boundary_ids (pt_neighbor, ns); 00083 00084 bool on_relevant_boundary = false; 00085 for (std::set<boundary_id_type>::const_iterator i = 00086 boundary_ids.begin(); i != boundary_ids.end(); ++i) 00087 if (std::find(bc_ids.begin(), bc_ids.end(), *i) 00088 != bc_ids.end()) 00089 on_relevant_boundary = true; 00090 00091 if (!on_relevant_boundary) 00092 continue; 00093 00094 if (!pt_neighbor->build_side(ns)->contains_point(p)) 00095 continue; 00096 00097 vertex_on_periodic_side = true; 00098 break; 00099 } 00100 00101 if (vertex_on_periodic_side) 00102 primary = pt_neighbor; 00103 } 00104 00105 return primary; 00106 } 00107 00108 // Find the "primary" element around a boundary edge: 00109 const Elem* primary_boundary_edge_neighbor 00110 (const Elem* elem, 00111 const Point& p1, 00112 const Point& p2, 00113 const BoundaryInfo& boundary_info, 00114 const std::set<boundary_id_type>& boundary_ids) 00115 { 00116 // If we don't find a better alternative, the user will have 00117 // provided the primary element 00118 const Elem *primary = elem; 00119 00120 std::set<const Elem*> edge_neighbors; 00121 elem->find_edge_neighbors(p1, p2, edge_neighbors); 00122 for (std::set<const Elem*>::const_iterator edge_neighbors_iter = 00123 edge_neighbors.begin(); 00124 edge_neighbors_iter != edge_neighbors.end(); ++edge_neighbors_iter) 00125 { 00126 const Elem* e_neighbor = *edge_neighbors_iter; 00127 00128 // If this edge neighbor isn't at least 00129 // as coarse as the current primary elem, or if it is at 00130 // the same level but has a lower id, then 00131 // we won't defer to it. 00132 if ((e_neighbor->level() > primary->level()) || 00133 (e_neighbor->level() == primary->level() && 00134 e_neighbor->id() < primary->id())) 00135 continue; 00136 00137 // Otherwise, we will defer to the edge neighbor, but only if 00138 // one of its sides is on this periodic boundary and that 00139 // side contains this edge 00140 bool vertex_on_periodic_side = false; 00141 for (unsigned int ns = 0; 00142 ns != e_neighbor->n_sides(); ++ns) 00143 { 00144 const std::vector<boundary_id_type>& bc_ids = 00145 boundary_info.boundary_ids (e_neighbor, ns); 00146 00147 bool on_relevant_boundary = false; 00148 for (std::set<boundary_id_type>::const_iterator i = 00149 boundary_ids.begin(); i != boundary_ids.end(); ++i) 00150 if (std::find(bc_ids.begin(), bc_ids.end(), *i) 00151 != bc_ids.end()) 00152 on_relevant_boundary = true; 00153 00154 if (!on_relevant_boundary) 00155 continue; 00156 00157 AutoPtr<Elem> periodic_side = e_neighbor->build_side(ns); 00158 if (!(periodic_side->contains_point(p1) && 00159 periodic_side->contains_point(p2))) 00160 continue; 00161 00162 vertex_on_periodic_side = true; 00163 break; 00164 } 00165 00166 if (vertex_on_periodic_side) 00167 primary = e_neighbor; 00168 } 00169 00170 return primary; 00171 } 00172 00173 } 00174 00175 namespace libMesh 00176 { 00177 00178 00179 00180 // ------------------------------------------------------------ 00181 // FEBase class members 00182 template <> 00183 AutoPtr<FEGenericBase<Real> > 00184 FEGenericBase<Real>::build (const unsigned int dim, 00185 const FEType& fet) 00186 { 00187 // The stupid AutoPtr<FEBase> ap(); return ap; 00188 // construct is required to satisfy IBM's xlC 00189 00190 switch (dim) 00191 { 00192 // 0D 00193 case 0: 00194 { 00195 switch (fet.family) 00196 { 00197 case CLOUGH: 00198 { 00199 AutoPtr<FEBase> ap(new FE<0,CLOUGH>(fet)); 00200 return ap; 00201 } 00202 00203 case HERMITE: 00204 { 00205 AutoPtr<FEBase> ap(new FE<0,HERMITE>(fet)); 00206 return ap; 00207 } 00208 00209 case LAGRANGE: 00210 { 00211 AutoPtr<FEBase> ap(new FE<0,LAGRANGE>(fet)); 00212 return ap; 00213 } 00214 00215 case L2_LAGRANGE: 00216 { 00217 AutoPtr<FEBase> ap(new FE<0,L2_LAGRANGE>(fet)); 00218 return ap; 00219 } 00220 00221 case HIERARCHIC: 00222 { 00223 AutoPtr<FEBase> ap(new FE<0,HIERARCHIC>(fet)); 00224 return ap; 00225 } 00226 00227 case L2_HIERARCHIC: 00228 { 00229 AutoPtr<FEBase> ap(new FE<0,L2_HIERARCHIC>(fet)); 00230 return ap; 00231 } 00232 00233 case MONOMIAL: 00234 { 00235 AutoPtr<FEBase> ap(new FE<0,MONOMIAL>(fet)); 00236 return ap; 00237 } 00238 00239 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00240 case SZABAB: 00241 { 00242 AutoPtr<FEBase> ap(new FE<0,SZABAB>(fet)); 00243 return ap; 00244 } 00245 00246 case BERNSTEIN: 00247 { 00248 AutoPtr<FEBase> ap(new FE<0,BERNSTEIN>(fet)); 00249 return ap; 00250 } 00251 #endif 00252 00253 case XYZ: 00254 { 00255 AutoPtr<FEBase> ap(new FEXYZ<0>(fet)); 00256 return ap; 00257 } 00258 00259 case SCALAR: 00260 { 00261 AutoPtr<FEBase> ap(new FEScalar<0>(fet)); 00262 return ap; 00263 } 00264 00265 default: 00266 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00267 libmesh_error(); 00268 } 00269 } 00270 // 1D 00271 case 1: 00272 { 00273 switch (fet.family) 00274 { 00275 case CLOUGH: 00276 { 00277 AutoPtr<FEBase> ap(new FE<1,CLOUGH>(fet)); 00278 return ap; 00279 } 00280 00281 case HERMITE: 00282 { 00283 AutoPtr<FEBase> ap(new FE<1,HERMITE>(fet)); 00284 return ap; 00285 } 00286 00287 case LAGRANGE: 00288 { 00289 AutoPtr<FEBase> ap(new FE<1,LAGRANGE>(fet)); 00290 return ap; 00291 } 00292 00293 case L2_LAGRANGE: 00294 { 00295 AutoPtr<FEBase> ap(new FE<1,L2_LAGRANGE>(fet)); 00296 return ap; 00297 } 00298 00299 case HIERARCHIC: 00300 { 00301 AutoPtr<FEBase> ap(new FE<1,HIERARCHIC>(fet)); 00302 return ap; 00303 } 00304 00305 case L2_HIERARCHIC: 00306 { 00307 AutoPtr<FEBase> ap(new FE<1,L2_HIERARCHIC>(fet)); 00308 return ap; 00309 } 00310 00311 case MONOMIAL: 00312 { 00313 AutoPtr<FEBase> ap(new FE<1,MONOMIAL>(fet)); 00314 return ap; 00315 } 00316 00317 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00318 case SZABAB: 00319 { 00320 AutoPtr<FEBase> ap(new FE<1,SZABAB>(fet)); 00321 return ap; 00322 } 00323 00324 case BERNSTEIN: 00325 { 00326 AutoPtr<FEBase> ap(new FE<1,BERNSTEIN>(fet)); 00327 return ap; 00328 } 00329 #endif 00330 00331 case XYZ: 00332 { 00333 AutoPtr<FEBase> ap(new FEXYZ<1>(fet)); 00334 return ap; 00335 } 00336 00337 case SCALAR: 00338 { 00339 AutoPtr<FEBase> ap(new FEScalar<1>(fet)); 00340 return ap; 00341 } 00342 00343 default: 00344 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00345 libmesh_error(); 00346 } 00347 } 00348 00349 00350 // 2D 00351 case 2: 00352 { 00353 switch (fet.family) 00354 { 00355 case CLOUGH: 00356 { 00357 AutoPtr<FEBase> ap(new FE<2,CLOUGH>(fet)); 00358 return ap; 00359 } 00360 00361 case HERMITE: 00362 { 00363 AutoPtr<FEBase> ap(new FE<2,HERMITE>(fet)); 00364 return ap; 00365 } 00366 00367 case LAGRANGE: 00368 { 00369 AutoPtr<FEBase> ap(new FE<2,LAGRANGE>(fet)); 00370 return ap; 00371 } 00372 00373 case L2_LAGRANGE: 00374 { 00375 AutoPtr<FEBase> ap(new FE<2,L2_LAGRANGE>(fet)); 00376 return ap; 00377 } 00378 00379 case HIERARCHIC: 00380 { 00381 AutoPtr<FEBase> ap(new FE<2,HIERARCHIC>(fet)); 00382 return ap; 00383 } 00384 00385 case L2_HIERARCHIC: 00386 { 00387 AutoPtr<FEBase> ap(new FE<2,L2_HIERARCHIC>(fet)); 00388 return ap; 00389 } 00390 00391 case MONOMIAL: 00392 { 00393 AutoPtr<FEBase> ap(new FE<2,MONOMIAL>(fet)); 00394 return ap; 00395 } 00396 00397 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00398 case SZABAB: 00399 { 00400 AutoPtr<FEBase> ap(new FE<2,SZABAB>(fet)); 00401 return ap; 00402 } 00403 00404 case BERNSTEIN: 00405 { 00406 AutoPtr<FEBase> ap(new FE<2,BERNSTEIN>(fet)); 00407 return ap; 00408 } 00409 #endif 00410 00411 case XYZ: 00412 { 00413 AutoPtr<FEBase> ap(new FEXYZ<2>(fet)); 00414 return ap; 00415 } 00416 00417 case SCALAR: 00418 { 00419 AutoPtr<FEBase> ap(new FEScalar<2>(fet)); 00420 return ap; 00421 } 00422 default: 00423 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00424 libmesh_error(); 00425 } 00426 } 00427 00428 00429 // 3D 00430 case 3: 00431 { 00432 switch (fet.family) 00433 { 00434 case CLOUGH: 00435 { 00436 libMesh::out << "ERROR: Clough-Tocher elements currently only support 1D and 2D" 00437 << std::endl; 00438 libmesh_error(); 00439 } 00440 00441 case HERMITE: 00442 { 00443 AutoPtr<FEBase> ap(new FE<3,HERMITE>(fet)); 00444 return ap; 00445 } 00446 00447 case LAGRANGE: 00448 { 00449 AutoPtr<FEBase> ap(new FE<3,LAGRANGE>(fet)); 00450 return ap; 00451 } 00452 00453 case L2_LAGRANGE: 00454 { 00455 AutoPtr<FEBase> ap(new FE<3,L2_LAGRANGE>(fet)); 00456 return ap; 00457 } 00458 00459 case HIERARCHIC: 00460 { 00461 AutoPtr<FEBase> ap(new FE<3,HIERARCHIC>(fet)); 00462 return ap; 00463 } 00464 00465 case L2_HIERARCHIC: 00466 { 00467 AutoPtr<FEBase> ap(new FE<3,L2_HIERARCHIC>(fet)); 00468 return ap; 00469 } 00470 00471 case MONOMIAL: 00472 { 00473 AutoPtr<FEBase> ap(new FE<3,MONOMIAL>(fet)); 00474 return ap; 00475 } 00476 00477 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 00478 case SZABAB: 00479 { 00480 AutoPtr<FEBase> ap(new FE<3,SZABAB>(fet)); 00481 return ap; 00482 } 00483 00484 case BERNSTEIN: 00485 { 00486 AutoPtr<FEBase> ap(new FE<3,BERNSTEIN>(fet)); 00487 return ap; 00488 } 00489 #endif 00490 00491 case XYZ: 00492 { 00493 AutoPtr<FEBase> ap(new FEXYZ<3>(fet)); 00494 return ap; 00495 } 00496 00497 case SCALAR: 00498 { 00499 AutoPtr<FEBase> ap(new FEScalar<3>(fet)); 00500 return ap; 00501 } 00502 00503 default: 00504 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00505 libmesh_error(); 00506 } 00507 } 00508 00509 default: 00510 libmesh_error(); 00511 } 00512 00513 libmesh_error(); 00514 AutoPtr<FEBase> ap(NULL); 00515 return ap; 00516 } 00517 00518 00519 00520 template <> 00521 AutoPtr<FEGenericBase<RealGradient> > 00522 FEGenericBase<RealGradient>::build (const unsigned int dim, 00523 const FEType& fet) 00524 { 00525 // The stupid AutoPtr<FEBase> ap(); return ap; 00526 // construct is required to satisfy IBM's xlC 00527 00528 switch (dim) 00529 { 00530 // 0D 00531 case 0: 00532 { 00533 switch (fet.family) 00534 { 00535 case LAGRANGE_VEC: 00536 { 00537 AutoPtr<FEVectorBase> ap( new FELagrangeVec<0>(fet) ); 00538 return ap; 00539 } 00540 default: 00541 { 00542 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00543 libmesh_error(); 00544 } 00545 } 00546 } 00547 case 1: 00548 { 00549 switch (fet.family) 00550 { 00551 case LAGRANGE_VEC: 00552 { 00553 AutoPtr<FEVectorBase> ap( new FELagrangeVec<1>(fet) ); 00554 return ap; 00555 } 00556 default: 00557 { 00558 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00559 libmesh_error(); 00560 } 00561 } 00562 } 00563 case 2: 00564 { 00565 switch (fet.family) 00566 { 00567 case LAGRANGE_VEC: 00568 { 00569 AutoPtr<FEVectorBase> ap( new FELagrangeVec<2>(fet) ); 00570 return ap; 00571 } 00572 case NEDELEC_ONE: 00573 { 00574 AutoPtr<FEVectorBase> ap( new FENedelecOne<2>(fet) ); 00575 return ap; 00576 } 00577 default: 00578 { 00579 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00580 libmesh_error(); 00581 } 00582 } 00583 } 00584 case 3: 00585 { 00586 switch (fet.family) 00587 { 00588 case LAGRANGE_VEC: 00589 { 00590 AutoPtr<FEVectorBase> ap( new FELagrangeVec<3>(fet) ); 00591 return ap; 00592 } 00593 case NEDELEC_ONE: 00594 { 00595 AutoPtr<FEVectorBase> ap( new FENedelecOne<3>(fet) ); 00596 return ap; 00597 } 00598 default: 00599 { 00600 libMesh::out << "ERROR: Bad FEType.family= " << fet.family << std::endl; 00601 libmesh_error(); 00602 } 00603 } 00604 } 00605 00606 default: 00607 libmesh_error(); 00608 00609 } // switch(dim) 00610 00611 libmesh_error(); 00612 AutoPtr<FEVectorBase> ap(NULL); 00613 return ap; 00614 } 00615 00616 00617 00618 00619 00620 00621 00622 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00623 00624 00625 template <> 00626 AutoPtr<FEGenericBase<Real> > 00627 FEGenericBase<Real>::build_InfFE (const unsigned int dim, 00628 const FEType& fet) 00629 { 00630 // The stupid AutoPtr<FEBase> ap(); return ap; 00631 // construct is required to satisfy IBM's xlC 00632 00633 switch (dim) 00634 { 00635 00636 // 1D 00637 case 1: 00638 { 00639 switch (fet.radial_family) 00640 { 00641 case INFINITE_MAP: 00642 { 00643 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00644 << " with FEFamily = " << fet.radial_family << std::endl; 00645 libmesh_error(); 00646 } 00647 00648 case JACOBI_20_00: 00649 { 00650 switch (fet.inf_map) 00651 { 00652 case CARTESIAN: 00653 { 00654 AutoPtr<FEBase> ap(new InfFE<1,JACOBI_20_00,CARTESIAN>(fet)); 00655 return ap; 00656 } 00657 default: 00658 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00659 << " with InfMapType = " << fet.inf_map << std::endl; 00660 libmesh_error(); 00661 } 00662 } 00663 00664 case JACOBI_30_00: 00665 { 00666 switch (fet.inf_map) 00667 { 00668 case CARTESIAN: 00669 { 00670 AutoPtr<FEBase> ap(new InfFE<1,JACOBI_30_00,CARTESIAN>(fet)); 00671 return ap; 00672 } 00673 default: 00674 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00675 << " with InfMapType = " << fet.inf_map << std::endl; 00676 libmesh_error(); 00677 } 00678 } 00679 00680 case LEGENDRE: 00681 { 00682 switch (fet.inf_map) 00683 { 00684 case CARTESIAN: 00685 { 00686 AutoPtr<FEBase> ap(new InfFE<1,LEGENDRE,CARTESIAN>(fet)); 00687 return ap; 00688 } 00689 default: 00690 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00691 << " with InfMapType = " << fet.inf_map << std::endl; 00692 libmesh_error(); 00693 } 00694 } 00695 00696 case LAGRANGE: 00697 { 00698 switch (fet.inf_map) 00699 { 00700 case CARTESIAN: 00701 { 00702 AutoPtr<FEBase> ap(new InfFE<1,LAGRANGE,CARTESIAN>(fet)); 00703 return ap; 00704 } 00705 default: 00706 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00707 << " with InfMapType = " << fet.inf_map << std::endl; 00708 libmesh_error(); 00709 } 00710 } 00711 00712 00713 00714 default: 00715 libMesh::err << "ERROR: Bad FEType.radial_family= " << fet.radial_family << std::endl; 00716 libmesh_error(); 00717 } 00718 00719 } 00720 00721 00722 00723 00724 // 2D 00725 case 2: 00726 { 00727 switch (fet.radial_family) 00728 { 00729 case INFINITE_MAP: 00730 { 00731 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00732 << " with FEFamily = " << fet.radial_family << std::endl; 00733 libmesh_error(); 00734 } 00735 00736 case JACOBI_20_00: 00737 { 00738 switch (fet.inf_map) 00739 { 00740 case CARTESIAN: 00741 { 00742 AutoPtr<FEBase> ap(new InfFE<2,JACOBI_20_00,CARTESIAN>(fet)); 00743 return ap; 00744 } 00745 default: 00746 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00747 << " with InfMapType = " << fet.inf_map << std::endl; 00748 libmesh_error(); 00749 } 00750 } 00751 00752 case JACOBI_30_00: 00753 { 00754 switch (fet.inf_map) 00755 { 00756 case CARTESIAN: 00757 { 00758 AutoPtr<FEBase> ap(new InfFE<2,JACOBI_30_00,CARTESIAN>(fet)); 00759 return ap; 00760 } 00761 default: 00762 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00763 << " with InfMapType = " << fet.inf_map << std::endl; 00764 libmesh_error(); 00765 } 00766 } 00767 00768 case LEGENDRE: 00769 { 00770 switch (fet.inf_map) 00771 { 00772 case CARTESIAN: 00773 { 00774 AutoPtr<FEBase> ap(new InfFE<2,LEGENDRE,CARTESIAN>(fet)); 00775 return ap; 00776 } 00777 default: 00778 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00779 << " with InfMapType = " << fet.inf_map << std::endl; 00780 libmesh_error(); 00781 } 00782 } 00783 00784 case LAGRANGE: 00785 { 00786 switch (fet.inf_map) 00787 { 00788 case CARTESIAN: 00789 { 00790 AutoPtr<FEBase> ap(new InfFE<2,LAGRANGE,CARTESIAN>(fet)); 00791 return ap; 00792 } 00793 default: 00794 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00795 << " with InfMapType = " << fet.inf_map << std::endl; 00796 libmesh_error(); 00797 } 00798 } 00799 00800 00801 00802 default: 00803 libMesh::err << "ERROR: Bad FEType.radial_family= " << fet.radial_family << std::endl; 00804 libmesh_error(); 00805 } 00806 00807 } 00808 00809 00810 00811 00812 // 3D 00813 case 3: 00814 { 00815 switch (fet.radial_family) 00816 { 00817 case INFINITE_MAP: 00818 { 00819 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00820 << " with FEFamily = " << fet.radial_family << std::endl; 00821 libmesh_error(); 00822 } 00823 00824 case JACOBI_20_00: 00825 { 00826 switch (fet.inf_map) 00827 { 00828 case CARTESIAN: 00829 { 00830 AutoPtr<FEBase> ap(new InfFE<3,JACOBI_20_00,CARTESIAN>(fet)); 00831 return ap; 00832 } 00833 default: 00834 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00835 << " with InfMapType = " << fet.inf_map << std::endl; 00836 libmesh_error(); 00837 } 00838 } 00839 00840 case JACOBI_30_00: 00841 { 00842 switch (fet.inf_map) 00843 { 00844 case CARTESIAN: 00845 { 00846 AutoPtr<FEBase> ap(new InfFE<3,JACOBI_30_00,CARTESIAN>(fet)); 00847 return ap; 00848 } 00849 default: 00850 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00851 << " with InfMapType = " << fet.inf_map << std::endl; 00852 libmesh_error(); 00853 } 00854 } 00855 00856 case LEGENDRE: 00857 { 00858 switch (fet.inf_map) 00859 { 00860 case CARTESIAN: 00861 { 00862 AutoPtr<FEBase> ap(new InfFE<3,LEGENDRE,CARTESIAN>(fet)); 00863 return ap; 00864 } 00865 default: 00866 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00867 << " with InfMapType = " << fet.inf_map << std::endl; 00868 libmesh_error(); 00869 } 00870 } 00871 00872 case LAGRANGE: 00873 { 00874 switch (fet.inf_map) 00875 { 00876 case CARTESIAN: 00877 { 00878 AutoPtr<FEBase> ap(new InfFE<3,LAGRANGE,CARTESIAN>(fet)); 00879 return ap; 00880 } 00881 default: 00882 libMesh::err << "ERROR: Don't build an infinite element " << std::endl 00883 << " with InfMapType = " << fet.inf_map << std::endl; 00884 libmesh_error(); 00885 } 00886 } 00887 00888 00889 00890 default: 00891 libMesh::err << "ERROR: Bad FEType.radial_family= " << fet.radial_family << std::endl; 00892 libmesh_error(); 00893 } 00894 } 00895 00896 default: 00897 libmesh_error(); 00898 } 00899 00900 libmesh_error(); 00901 AutoPtr<FEBase> ap(NULL); 00902 return ap; 00903 } 00904 00905 00906 00907 template <> 00908 AutoPtr<FEGenericBase<RealGradient> > 00909 FEGenericBase<RealGradient>::build_InfFE (const unsigned int, 00910 const FEType& ) 00911 { 00912 // No vector types defined... YET. 00913 libmesh_error(); 00914 AutoPtr<FEVectorBase> ap(NULL); 00915 return ap; 00916 } 00917 00918 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00919 00920 00921 template <typename OutputType> 00922 void FEGenericBase<OutputType> ::compute_shape_functions (const Elem* elem, 00923 const std::vector<Point>& qp) 00924 { 00925 //------------------------------------------------------------------------- 00926 // Compute the shape function values (and derivatives) 00927 // at the Quadrature points. Note that the actual values 00928 // have already been computed via init_shape_functions 00929 00930 // Start logging the shape function computation 00931 START_LOG("compute_shape_functions()", "FE"); 00932 00933 calculations_started = true; 00934 00935 // If the user forgot to request anything, we'll be safe and 00936 // calculate everything: 00937 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00938 if (!calculate_phi && !calculate_dphi && !calculate_d2phi && !calculate_curl_phi && !calculate_div_phi) 00939 { 00940 calculate_phi = calculate_dphi = calculate_d2phi = true; 00941 // Only compute curl, div for vector-valued elements 00942 if( TypesEqual<OutputType,RealGradient>::value ) 00943 { 00944 calculate_curl_phi = true; 00945 calculate_div_phi = true; 00946 } 00947 } 00948 #else 00949 if (!calculate_phi && !calculate_dphi && !calculate_curl_phi && !calculate_div_phi) 00950 { 00951 calculate_phi = calculate_dphi = true; 00952 // Only compute curl for vector-valued elements 00953 if( TypesEqual<OutputType,RealGradient>::value ) 00954 { 00955 calculate_curl_phi = true; 00956 calculate_div_phi = true; 00957 } 00958 } 00959 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 00960 00961 00962 if( calculate_phi ) 00963 this->_fe_trans->map_phi( this->dim, elem, qp, (*this), this->phi ); 00964 00965 if( calculate_dphi ) 00966 this->_fe_trans->map_dphi( this->dim, elem, qp, (*this), this->dphi, 00967 this->dphidx, this->dphidy, this->dphidz ); 00968 00969 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00970 if( calculate_d2phi ) 00971 this->_fe_trans->map_d2phi( this->dim, elem, qp, (*this), this->d2phi, 00972 this->d2phidx2, this->d2phidxdy, this->d2phidxdz, 00973 this->d2phidy2, this->d2phidydz, this->d2phidz2 ); 00974 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES 00975 00976 // Only compute curl for vector-valued elements 00977 if( calculate_curl_phi && TypesEqual<OutputType,RealGradient>::value ) 00978 this->_fe_trans->map_curl( this->dim, elem, qp, (*this), this->curl_phi ); 00979 00980 // Only compute div for vector-valued elements 00981 if( calculate_div_phi && TypesEqual<OutputType,RealGradient>::value ) 00982 this->_fe_trans->map_div( this->dim, elem, qp, (*this), this->div_phi ); 00983 00984 // Stop logging the shape function computation 00985 STOP_LOG("compute_shape_functions()", "FE"); 00986 } 00987 00988 00989 template <typename OutputType> 00990 void FEGenericBase<OutputType>::print_phi(std::ostream& os) const 00991 { 00992 for (unsigned int i=0; i<phi.size(); ++i) 00993 for (unsigned int j=0; j<phi[i].size(); ++j) 00994 os << " phi[" << i << "][" << j << "]=" << phi[i][j] << std::endl; 00995 } 00996 00997 00998 00999 01000 template <typename OutputType> 01001 void FEGenericBase<OutputType>::print_dphi(std::ostream& os) const 01002 { 01003 for (unsigned int i=0; i<dphi.size(); ++i) 01004 for (unsigned int j=0; j<dphi[i].size(); ++j) 01005 os << " dphi[" << i << "][" << j << "]=" << dphi[i][j]; 01006 } 01007 01008 01009 01010 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 01011 01012 01013 template <typename OutputType> 01014 void FEGenericBase<OutputType>::print_d2phi(std::ostream& os) const 01015 { 01016 for (unsigned int i=0; i<dphi.size(); ++i) 01017 for (unsigned int j=0; j<dphi[i].size(); ++j) 01018 os << " d2phi[" << i << "][" << j << "]=" << d2phi[i][j]; 01019 } 01020 01021 #endif 01022 01023 01024 01025 #ifdef LIBMESH_ENABLE_AMR 01026 01027 template <typename OutputType> 01028 void 01029 FEGenericBase<OutputType>::coarsened_dof_values(const NumericVector<Number> &old_vector, 01030 const DofMap &dof_map, 01031 const Elem *elem, 01032 DenseVector<Number> &Ue, 01033 const unsigned int var, 01034 const bool use_old_dof_indices) 01035 { 01036 // Side/edge local DOF indices 01037 std::vector<unsigned int> new_side_dofs, old_side_dofs; 01038 01039 // FIXME: what about 2D shells in 3D space? 01040 unsigned int dim = elem->dim(); 01041 01042 // We use local FE objects for now 01043 // FIXME: we should use more, external objects instead for efficiency 01044 const FEType& base_fe_type = dof_map.variable_type(var); 01045 AutoPtr<FEGenericBase<OutputShape> > fe 01046 (FEGenericBase<OutputShape>::build(dim, base_fe_type)); 01047 AutoPtr<FEGenericBase<OutputShape> > fe_coarse 01048 (FEGenericBase<OutputShape>::build(dim, base_fe_type)); 01049 01050 AutoPtr<QBase> qrule (base_fe_type.default_quadrature_rule(dim)); 01051 AutoPtr<QBase> qedgerule (base_fe_type.default_quadrature_rule(1)); 01052 AutoPtr<QBase> qsiderule (base_fe_type.default_quadrature_rule(dim-1)); 01053 std::vector<Point> coarse_qpoints; 01054 01055 // The values of the shape functions at the quadrature 01056 // points 01057 const std::vector<std::vector<OutputShape> >& phi_values = 01058 fe->get_phi(); 01059 const std::vector<std::vector<OutputShape> >& phi_coarse = 01060 fe_coarse->get_phi(); 01061 01062 // The gradients of the shape functions at the quadrature 01063 // points on the child element. 01064 const std::vector<std::vector<OutputGradient> > *dphi_values = 01065 NULL; 01066 const std::vector<std::vector<OutputGradient> > *dphi_coarse = 01067 NULL; 01068 01069 const FEContinuity cont = fe->get_continuity(); 01070 01071 if (cont == C_ONE) 01072 { 01073 const std::vector<std::vector<OutputGradient> >& 01074 ref_dphi_values = fe->get_dphi(); 01075 dphi_values = &ref_dphi_values; 01076 const std::vector<std::vector<OutputGradient> >& 01077 ref_dphi_coarse = fe_coarse->get_dphi(); 01078 dphi_coarse = &ref_dphi_coarse; 01079 } 01080 01081 // The Jacobian * quadrature weight at the quadrature points 01082 const std::vector<Real>& JxW = 01083 fe->get_JxW(); 01084 01085 // The XYZ locations of the quadrature points on the 01086 // child element 01087 const std::vector<Point>& xyz_values = 01088 fe->get_xyz(); 01089 01090 01091 01092 FEType fe_type = base_fe_type, temp_fe_type; 01093 const ElemType elem_type = elem->type(); 01094 fe_type.order = static_cast<Order>(fe_type.order + 01095 elem->max_descendant_p_level()); 01096 01097 // Number of nodes on parent element 01098 const unsigned int n_nodes = elem->n_nodes(); 01099 01100 // Number of dofs on parent element 01101 const unsigned int new_n_dofs = 01102 FEInterface::n_dofs(dim, fe_type, elem_type); 01103 01104 // Fixed vs. free DoFs on edge/face projections 01105 std::vector<char> dof_is_fixed(new_n_dofs, false); // bools 01106 std::vector<int> free_dof(new_n_dofs, 0); 01107 01108 DenseMatrix<Real> Ke; 01109 DenseVector<Number> Fe; 01110 Ue.resize(new_n_dofs); Ue.zero(); 01111 01112 01113 // When coarsening, in general, we need a series of 01114 // projections to ensure a unique and continuous 01115 // solution. We start by interpolating nodes, then 01116 // hold those fixed and project edges, then 01117 // hold those fixed and project faces, then 01118 // hold those fixed and project interiors 01119 01120 // Copy node values first 01121 { 01122 std::vector<dof_id_type> node_dof_indices; 01123 if (use_old_dof_indices) 01124 dof_map.old_dof_indices (elem, node_dof_indices, var); 01125 else 01126 dof_map.dof_indices (elem, node_dof_indices, var); 01127 01128 unsigned int current_dof = 0; 01129 for (unsigned int n=0; n!= n_nodes; ++n) 01130 { 01131 // FIXME: this should go through the DofMap, 01132 // not duplicate dof_indices code badly! 01133 const unsigned int my_nc = 01134 FEInterface::n_dofs_at_node (dim, fe_type, 01135 elem_type, n); 01136 if (!elem->is_vertex(n)) 01137 { 01138 current_dof += my_nc; 01139 continue; 01140 } 01141 01142 temp_fe_type = base_fe_type; 01143 // We're assuming here that child n shares vertex n, 01144 // which is wrong on non-simplices right now 01145 // ... but this code isn't necessary except on elements 01146 // where p refinement creates more vertex dofs; we have 01147 // no such elements yet. 01148 /* 01149 if (elem->child(n)->p_level() < elem->p_level()) 01150 { 01151 temp_fe_type.order = 01152 static_cast<Order>(temp_fe_type.order + 01153 elem->child(n)->p_level()); 01154 } 01155 */ 01156 const unsigned int nc = 01157 FEInterface::n_dofs_at_node (dim, temp_fe_type, 01158 elem_type, n); 01159 for (unsigned int i=0; i!= nc; ++i) 01160 { 01161 Ue(current_dof) = 01162 old_vector(node_dof_indices[current_dof]); 01163 dof_is_fixed[current_dof] = true; 01164 current_dof++; 01165 } 01166 } 01167 } 01168 01169 // In 3D, project any edge values next 01170 if (dim > 2 && cont != DISCONTINUOUS) 01171 for (unsigned int e=0; e != elem->n_edges(); ++e) 01172 { 01173 FEInterface::dofs_on_edge(elem, dim, fe_type, 01174 e, new_side_dofs); 01175 01176 // Some edge dofs are on nodes and already 01177 // fixed, others are free to calculate 01178 unsigned int free_dofs = 0; 01179 for (unsigned int i=0; i != 01180 new_side_dofs.size(); ++i) 01181 if (!dof_is_fixed[new_side_dofs[i]]) 01182 free_dof[free_dofs++] = i; 01183 Ke.resize (free_dofs, free_dofs); Ke.zero(); 01184 Fe.resize (free_dofs); Fe.zero(); 01185 // The new edge coefficients 01186 DenseVector<Number> Uedge(free_dofs); 01187 01188 // Add projection terms from each child sharing 01189 // this edge 01190 for (unsigned int c=0; c != elem->n_children(); 01191 ++c) 01192 { 01193 if (!elem->is_child_on_edge(c,e)) 01194 continue; 01195 Elem *child = elem->child(c); 01196 01197 std::vector<dof_id_type> child_dof_indices; 01198 if (use_old_dof_indices) 01199 dof_map.old_dof_indices (child, 01200 child_dof_indices, var); 01201 else 01202 dof_map.dof_indices (child, 01203 child_dof_indices, var); 01204 const unsigned int child_n_dofs = 01205 libmesh_cast_int<unsigned int> 01206 (child_dof_indices.size()); 01207 01208 temp_fe_type = base_fe_type; 01209 temp_fe_type.order = 01210 static_cast<Order>(temp_fe_type.order + 01211 child->p_level()); 01212 01213 FEInterface::dofs_on_edge(child, dim, 01214 temp_fe_type, e, old_side_dofs); 01215 01216 // Initialize both child and parent FE data 01217 // on the child's edge 01218 fe->attach_quadrature_rule (qedgerule.get()); 01219 fe->edge_reinit (child, e); 01220 const unsigned int n_qp = qedgerule->n_points(); 01221 01222 FEInterface::inverse_map (dim, fe_type, elem, 01223 xyz_values, coarse_qpoints); 01224 01225 fe_coarse->reinit(elem, &coarse_qpoints); 01226 01227 // Loop over the quadrature points 01228 for (unsigned int qp=0; qp<n_qp; qp++) 01229 { 01230 // solution value at the quadrature point 01231 OutputNumber fineval = libMesh::zero; 01232 // solution grad at the quadrature point 01233 OutputNumberGradient finegrad; 01234 01235 // Sum the solution values * the DOF 01236 // values at the quadrature point to 01237 // get the solution value and gradient. 01238 for (unsigned int i=0; i<child_n_dofs; 01239 i++) 01240 { 01241 fineval += 01242 (old_vector(child_dof_indices[i])* 01243 phi_values[i][qp]); 01244 if (cont == C_ONE) 01245 finegrad += (*dphi_values)[i][qp] * 01246 old_vector(child_dof_indices[i]); 01247 } 01248 01249 // Form edge projection matrix 01250 for (unsigned int sidei=0, freei=0; 01251 sidei != new_side_dofs.size(); 01252 ++sidei) 01253 { 01254 unsigned int i = new_side_dofs[sidei]; 01255 // fixed DoFs aren't test functions 01256 if (dof_is_fixed[i]) 01257 continue; 01258 for (unsigned int sidej=0, freej=0; 01259 sidej != new_side_dofs.size(); 01260 ++sidej) 01261 { 01262 unsigned int j = 01263 new_side_dofs[sidej]; 01264 if (dof_is_fixed[j]) 01265 Fe(freei) -= 01266 TensorTools::inner_product(phi_coarse[i][qp], 01267 phi_coarse[j][qp]) * 01268 JxW[qp] * Ue(j); 01269 else 01270 Ke(freei,freej) += 01271 TensorTools::inner_product(phi_coarse[i][qp], 01272 phi_coarse[j][qp]) * 01273 JxW[qp]; 01274 if (cont == C_ONE) 01275 { 01276 if (dof_is_fixed[j]) 01277 Fe(freei) -= 01278 TensorTools::inner_product((*dphi_coarse)[i][qp], 01279 (*dphi_coarse)[j][qp]) * 01280 JxW[qp] * Ue(j); 01281 else 01282 Ke(freei,freej) += 01283 TensorTools::inner_product((*dphi_coarse)[i][qp], 01284 (*dphi_coarse)[j][qp]) * 01285 JxW[qp]; 01286 } 01287 if (!dof_is_fixed[j]) 01288 freej++; 01289 } 01290 Fe(freei) += TensorTools::inner_product(phi_coarse[i][qp], 01291 fineval) * JxW[qp]; 01292 if (cont == C_ONE) 01293 Fe(freei) += 01294 TensorTools::inner_product(finegrad, (*dphi_coarse)[i][qp]) * JxW[qp]; 01295 freei++; 01296 } 01297 } 01298 } 01299 Ke.cholesky_solve(Fe, Uedge); 01300 01301 // Transfer new edge solutions to element 01302 for (unsigned int i=0; i != free_dofs; ++i) 01303 { 01304 Number &ui = Ue(new_side_dofs[free_dof[i]]); 01305 libmesh_assert(std::abs(ui) < TOLERANCE || 01306 std::abs(ui - Uedge(i)) < TOLERANCE); 01307 ui = Uedge(i); 01308 dof_is_fixed[new_side_dofs[free_dof[i]]] = 01309 true; 01310 } 01311 } 01312 01313 // Project any side values (edges in 2D, faces in 3D) 01314 if (dim > 1 && cont != DISCONTINUOUS) 01315 for (unsigned int s=0; s != elem->n_sides(); ++s) 01316 { 01317 FEInterface::dofs_on_side(elem, dim, fe_type, 01318 s, new_side_dofs); 01319 01320 // Some side dofs are on nodes/edges and already 01321 // fixed, others are free to calculate 01322 unsigned int free_dofs = 0; 01323 for (unsigned int i=0; i != 01324 new_side_dofs.size(); ++i) 01325 if (!dof_is_fixed[new_side_dofs[i]]) 01326 free_dof[free_dofs++] = i; 01327 Ke.resize (free_dofs, free_dofs); Ke.zero(); 01328 Fe.resize (free_dofs); Fe.zero(); 01329 // The new side coefficients 01330 DenseVector<Number> Uside(free_dofs); 01331 01332 // Add projection terms from each child sharing 01333 // this side 01334 for (unsigned int c=0; c != elem->n_children(); 01335 ++c) 01336 { 01337 if (!elem->is_child_on_side(c,s)) 01338 continue; 01339 Elem *child = elem->child(c); 01340 01341 std::vector<dof_id_type> child_dof_indices; 01342 if (use_old_dof_indices) 01343 dof_map.old_dof_indices (child, 01344 child_dof_indices, var); 01345 else 01346 dof_map.dof_indices (child, 01347 child_dof_indices, var); 01348 const unsigned int child_n_dofs = 01349 libmesh_cast_int<unsigned int> 01350 (child_dof_indices.size()); 01351 01352 temp_fe_type = base_fe_type; 01353 temp_fe_type.order = 01354 static_cast<Order>(temp_fe_type.order + 01355 child->p_level()); 01356 01357 FEInterface::dofs_on_side(child, dim, 01358 temp_fe_type, s, old_side_dofs); 01359 01360 // Initialize both child and parent FE data 01361 // on the child's side 01362 fe->attach_quadrature_rule (qsiderule.get()); 01363 fe->reinit (child, s); 01364 const unsigned int n_qp = qsiderule->n_points(); 01365 01366 FEInterface::inverse_map (dim, fe_type, elem, 01367 xyz_values, coarse_qpoints); 01368 01369 fe_coarse->reinit(elem, &coarse_qpoints); 01370 01371 // Loop over the quadrature points 01372 for (unsigned int qp=0; qp<n_qp; qp++) 01373 { 01374 // solution value at the quadrature point 01375 OutputNumber fineval = libMesh::zero; 01376 // solution grad at the quadrature point 01377 OutputNumberGradient finegrad; 01378 01379 // Sum the solution values * the DOF 01380 // values at the quadrature point to 01381 // get the solution value and gradient. 01382 for (unsigned int i=0; i<child_n_dofs; 01383 i++) 01384 { 01385 fineval += 01386 old_vector(child_dof_indices[i]) * 01387 phi_values[i][qp]; 01388 if (cont == C_ONE) 01389 finegrad += (*dphi_values)[i][qp] * 01390 old_vector(child_dof_indices[i]); 01391 } 01392 01393 // Form side projection matrix 01394 for (unsigned int sidei=0, freei=0; 01395 sidei != new_side_dofs.size(); 01396 ++sidei) 01397 { 01398 unsigned int i = new_side_dofs[sidei]; 01399 // fixed DoFs aren't test functions 01400 if (dof_is_fixed[i]) 01401 continue; 01402 for (unsigned int sidej=0, freej=0; 01403 sidej != new_side_dofs.size(); 01404 ++sidej) 01405 { 01406 unsigned int j = 01407 new_side_dofs[sidej]; 01408 if (dof_is_fixed[j]) 01409 Fe(freei) -= 01410 TensorTools::inner_product(phi_coarse[i][qp], 01411 phi_coarse[j][qp]) * 01412 JxW[qp] * Ue(j); 01413 else 01414 Ke(freei,freej) += 01415 TensorTools::inner_product(phi_coarse[i][qp], 01416 phi_coarse[j][qp]) * 01417 JxW[qp]; 01418 if (cont == C_ONE) 01419 { 01420 if (dof_is_fixed[j]) 01421 Fe(freei) -= 01422 TensorTools::inner_product((*dphi_coarse)[i][qp], 01423 (*dphi_coarse)[j][qp]) * 01424 JxW[qp] * Ue(j); 01425 else 01426 Ke(freei,freej) += 01427 TensorTools::inner_product((*dphi_coarse)[i][qp], 01428 (*dphi_coarse)[j][qp]) * 01429 JxW[qp]; 01430 } 01431 if (!dof_is_fixed[j]) 01432 freej++; 01433 } 01434 Fe(freei) += TensorTools::inner_product(fineval, phi_coarse[i][qp]) * JxW[qp]; 01435 if (cont == C_ONE) 01436 Fe(freei) += 01437 TensorTools::inner_product(finegrad, (*dphi_coarse)[i][qp]) * JxW[qp]; 01438 freei++; 01439 } 01440 } 01441 } 01442 Ke.cholesky_solve(Fe, Uside); 01443 01444 // Transfer new side solutions to element 01445 for (unsigned int i=0; i != free_dofs; ++i) 01446 { 01447 Number &ui = Ue(new_side_dofs[free_dof[i]]); 01448 libmesh_assert(std::abs(ui) < TOLERANCE || 01449 std::abs(ui - Uside(i)) < TOLERANCE); 01450 ui = Uside(i); 01451 dof_is_fixed[new_side_dofs[free_dof[i]]] = 01452 true; 01453 } 01454 } 01455 01456 // Project the interior values, finally 01457 01458 // Some interior dofs are on nodes/edges/sides and 01459 // already fixed, others are free to calculate 01460 unsigned int free_dofs = 0; 01461 for (unsigned int i=0; i != new_n_dofs; ++i) 01462 if (!dof_is_fixed[i]) 01463 free_dof[free_dofs++] = i; 01464 Ke.resize (free_dofs, free_dofs); Ke.zero(); 01465 Fe.resize (free_dofs); Fe.zero(); 01466 // The new interior coefficients 01467 DenseVector<Number> Uint(free_dofs); 01468 01469 // Add projection terms from each child 01470 for (unsigned int c=0; c != elem->n_children(); ++c) 01471 { 01472 Elem *child = elem->child(c); 01473 01474 std::vector<dof_id_type> child_dof_indices; 01475 if (use_old_dof_indices) 01476 dof_map.old_dof_indices (child, 01477 child_dof_indices, var); 01478 else 01479 dof_map.dof_indices (child, 01480 child_dof_indices, var); 01481 const unsigned int child_n_dofs = 01482 libmesh_cast_int<unsigned int> 01483 (child_dof_indices.size()); 01484 01485 // Initialize both child and parent FE data 01486 // on the child's quadrature points 01487 fe->attach_quadrature_rule (qrule.get()); 01488 fe->reinit (child); 01489 const unsigned int n_qp = qrule->n_points(); 01490 01491 FEInterface::inverse_map (dim, fe_type, elem, 01492 xyz_values, coarse_qpoints); 01493 01494 fe_coarse->reinit(elem, &coarse_qpoints); 01495 01496 // Loop over the quadrature points 01497 for (unsigned int qp=0; qp<n_qp; qp++) 01498 { 01499 // solution value at the quadrature point 01500 OutputNumber fineval = libMesh::zero; 01501 // solution grad at the quadrature point 01502 OutputNumberGradient finegrad; 01503 01504 // Sum the solution values * the DOF 01505 // values at the quadrature point to 01506 // get the solution value and gradient. 01507 for (unsigned int i=0; i<child_n_dofs; i++) 01508 { 01509 fineval += 01510 (old_vector(child_dof_indices[i]) * 01511 phi_values[i][qp]); 01512 if (cont == C_ONE) 01513 finegrad += (*dphi_values)[i][qp] * 01514 old_vector(child_dof_indices[i]); 01515 } 01516 01517 // Form interior projection matrix 01518 for (unsigned int i=0, freei=0; 01519 i != new_n_dofs; ++i) 01520 { 01521 // fixed DoFs aren't test functions 01522 if (dof_is_fixed[i]) 01523 continue; 01524 for (unsigned int j=0, freej=0; j != 01525 new_n_dofs; ++j) 01526 { 01527 if (dof_is_fixed[j]) 01528 Fe(freei) -= 01529 TensorTools::inner_product(phi_coarse[i][qp], 01530 phi_coarse[j][qp]) * 01531 JxW[qp] * Ue(j); 01532 else 01533 Ke(freei,freej) += 01534 TensorTools::inner_product(phi_coarse[i][qp], 01535 phi_coarse[j][qp]) * 01536 JxW[qp]; 01537 if (cont == C_ONE) 01538 { 01539 if (dof_is_fixed[j]) 01540 Fe(freei) -= 01541 TensorTools::inner_product((*dphi_coarse)[i][qp], 01542 (*dphi_coarse)[j][qp]) * 01543 JxW[qp] * Ue(j); 01544 else 01545 Ke(freei,freej) += 01546 TensorTools::inner_product((*dphi_coarse)[i][qp], 01547 (*dphi_coarse)[j][qp]) * 01548 JxW[qp]; 01549 } 01550 if (!dof_is_fixed[j]) 01551 freej++; 01552 } 01553 Fe(freei) += TensorTools::inner_product(phi_coarse[i][qp], fineval) * 01554 JxW[qp]; 01555 if (cont == C_ONE) 01556 Fe(freei) += TensorTools::inner_product(finegrad, (*dphi_coarse)[i][qp]) * JxW[qp]; 01557 freei++; 01558 } 01559 } 01560 } 01561 Ke.cholesky_solve(Fe, Uint); 01562 01563 // Transfer new interior solutions to element 01564 for (unsigned int i=0; i != free_dofs; ++i) 01565 { 01566 Number &ui = Ue(free_dof[i]); 01567 libmesh_assert(std::abs(ui) < TOLERANCE || 01568 std::abs(ui - Uint(i)) < TOLERANCE); 01569 ui = Uint(i); 01570 dof_is_fixed[free_dof[i]] = true; 01571 } 01572 01573 // Make sure every DoF got reached! 01574 for (unsigned int i=0; i != new_n_dofs; ++i) 01575 libmesh_assert(dof_is_fixed[i]); 01576 } 01577 01578 01579 01580 template <typename OutputType> 01581 void 01582 FEGenericBase<OutputType>::compute_proj_constraints (DofConstraints &constraints, 01583 DofMap &dof_map, 01584 const unsigned int variable_number, 01585 const Elem* elem) 01586 { 01587 libmesh_assert(elem); 01588 01589 const unsigned int Dim = elem->dim(); 01590 01591 // Only constrain elements in 2,3D. 01592 if (Dim == 1) 01593 return; 01594 01595 // Only constrain active elements with this method 01596 if (!elem->active()) 01597 return; 01598 01599 const FEType& base_fe_type = dof_map.variable_type(variable_number); 01600 01601 // Construct FE objects for this element and its neighbors. 01602 AutoPtr<FEGenericBase<OutputShape> > my_fe 01603 (FEGenericBase<OutputShape>::build(Dim, base_fe_type)); 01604 const FEContinuity cont = my_fe->get_continuity(); 01605 01606 // We don't need to constrain discontinuous elements 01607 if (cont == DISCONTINUOUS) 01608 return; 01609 libmesh_assert (cont == C_ZERO || cont == C_ONE); 01610 01611 AutoPtr<FEGenericBase<OutputShape> > neigh_fe 01612 (FEGenericBase<OutputShape>::build(Dim, base_fe_type)); 01613 01614 QGauss my_qface(Dim-1, base_fe_type.default_quadrature_order()); 01615 my_fe->attach_quadrature_rule (&my_qface); 01616 std::vector<Point> neigh_qface; 01617 01618 const std::vector<Real>& JxW = my_fe->get_JxW(); 01619 const std::vector<Point>& q_point = my_fe->get_xyz(); 01620 const std::vector<std::vector<OutputShape> >& phi = my_fe->get_phi(); 01621 const std::vector<std::vector<OutputShape> >& neigh_phi = 01622 neigh_fe->get_phi(); 01623 const std::vector<Point> *face_normals = NULL; 01624 const std::vector<std::vector<OutputGradient> > *dphi = NULL; 01625 const std::vector<std::vector<OutputGradient> > *neigh_dphi = NULL; 01626 01627 std::vector<dof_id_type> my_dof_indices, neigh_dof_indices; 01628 std::vector<unsigned int> my_side_dofs, neigh_side_dofs; 01629 01630 if (cont != C_ZERO) 01631 { 01632 const std::vector<Point>& ref_face_normals = 01633 my_fe->get_normals(); 01634 face_normals = &ref_face_normals; 01635 const std::vector<std::vector<OutputGradient> >& ref_dphi = 01636 my_fe->get_dphi(); 01637 dphi = &ref_dphi; 01638 const std::vector<std::vector<OutputGradient> >& ref_neigh_dphi = 01639 neigh_fe->get_dphi(); 01640 neigh_dphi = &ref_neigh_dphi; 01641 } 01642 01643 DenseMatrix<Real> Ke; 01644 DenseVector<Real> Fe; 01645 std::vector<DenseVector<Real> > Ue; 01646 01647 // Look at the element faces. Check to see if we need to 01648 // build constraints. 01649 for (unsigned int s=0; s<elem->n_sides(); s++) 01650 if (elem->neighbor(s) != NULL) 01651 { 01652 // Get pointers to the element's neighbor. 01653 const Elem* neigh = elem->neighbor(s); 01654 01655 // h refinement constraints: 01656 // constrain dofs shared between 01657 // this element and ones coarser 01658 // than this element. 01659 if (neigh->level() < elem->level()) 01660 { 01661 unsigned int s_neigh = neigh->which_neighbor_am_i(elem); 01662 libmesh_assert_less (s_neigh, neigh->n_neighbors()); 01663 01664 // Find the minimum p level; we build the h constraint 01665 // matrix with this and then constrain away all higher p 01666 // DoFs. 01667 libmesh_assert(neigh->active()); 01668 const unsigned int min_p_level = 01669 std::min(elem->p_level(), neigh->p_level()); 01670 01671 // we may need to make the FE objects reinit with the 01672 // minimum shared p_level 01673 // FIXME - I hate using const_cast<> and avoiding 01674 // accessor functions; there's got to be a 01675 // better way to do this! 01676 const unsigned int old_elem_level = elem->p_level(); 01677 if (old_elem_level != min_p_level) 01678 (const_cast<Elem *>(elem))->hack_p_level(min_p_level); 01679 const unsigned int old_neigh_level = neigh->p_level(); 01680 if (old_neigh_level != min_p_level) 01681 (const_cast<Elem *>(neigh))->hack_p_level(min_p_level); 01682 01683 my_fe->reinit(elem, s); 01684 01685 // This function gets called element-by-element, so there 01686 // will be a lot of memory allocation going on. We can 01687 // at least minimize this for the case of the dof indices 01688 // by efficiently preallocating the requisite storage. 01689 // n_nodes is not necessarily n_dofs, but it is better 01690 // than nothing! 01691 my_dof_indices.reserve (elem->n_nodes()); 01692 neigh_dof_indices.reserve (neigh->n_nodes()); 01693 01694 dof_map.dof_indices (elem, my_dof_indices, 01695 variable_number); 01696 dof_map.dof_indices (neigh, neigh_dof_indices, 01697 variable_number); 01698 01699 const unsigned int n_qp = my_qface.n_points(); 01700 01701 FEInterface::inverse_map (Dim, base_fe_type, neigh, 01702 q_point, neigh_qface); 01703 01704 neigh_fe->reinit(neigh, &neigh_qface); 01705 01706 // We're only concerned with DOFs whose values (and/or first 01707 // derivatives for C1 elements) are supported on side nodes 01708 FEInterface::dofs_on_side(elem, Dim, base_fe_type, s, my_side_dofs); 01709 FEInterface::dofs_on_side(neigh, Dim, base_fe_type, s_neigh, neigh_side_dofs); 01710 01711 // We're done with functions that examine Elem::p_level(), 01712 // so let's unhack those levels 01713 if (elem->p_level() != old_elem_level) 01714 (const_cast<Elem *>(elem))->hack_p_level(old_elem_level); 01715 if (neigh->p_level() != old_neigh_level) 01716 (const_cast<Elem *>(neigh))->hack_p_level(old_neigh_level); 01717 01718 const unsigned int n_side_dofs = 01719 libmesh_cast_int<unsigned int>(my_side_dofs.size()); 01720 libmesh_assert_equal_to (n_side_dofs, neigh_side_dofs.size()); 01721 01722 Ke.resize (n_side_dofs, n_side_dofs); 01723 Ue.resize(n_side_dofs); 01724 01725 // Form the projection matrix, (inner product of fine basis 01726 // functions against fine test functions) 01727 for (unsigned int is = 0; is != n_side_dofs; ++is) 01728 { 01729 const unsigned int i = my_side_dofs[is]; 01730 for (unsigned int js = 0; js != n_side_dofs; ++js) 01731 { 01732 const unsigned int j = my_side_dofs[js]; 01733 for (unsigned int qp = 0; qp != n_qp; ++qp) 01734 { 01735 Ke(is,js) += JxW[qp] * TensorTools::inner_product(phi[i][qp], phi[j][qp]); 01736 if (cont != C_ZERO) 01737 Ke(is,js) += JxW[qp] * 01738 TensorTools::inner_product((*dphi)[i][qp] * 01739 (*face_normals)[qp], 01740 (*dphi)[j][qp] * 01741 (*face_normals)[qp]); 01742 } 01743 } 01744 } 01745 01746 // Form the right hand sides, (inner product of coarse basis 01747 // functions against fine test functions) 01748 for (unsigned int is = 0; is != n_side_dofs; ++is) 01749 { 01750 const unsigned int i = neigh_side_dofs[is]; 01751 Fe.resize (n_side_dofs); 01752 for (unsigned int js = 0; js != n_side_dofs; ++js) 01753 { 01754 const unsigned int j = my_side_dofs[js]; 01755 for (unsigned int qp = 0; qp != n_qp; ++qp) 01756 { 01757 Fe(js) += JxW[qp] * 01758 TensorTools::inner_product(neigh_phi[i][qp], 01759 phi[j][qp]); 01760 if (cont != C_ZERO) 01761 Fe(js) += JxW[qp] * 01762 TensorTools::inner_product((*neigh_dphi)[i][qp] * 01763 (*face_normals)[qp], 01764 (*dphi)[j][qp] * 01765 (*face_normals)[qp]); 01766 } 01767 } 01768 Ke.cholesky_solve(Fe, Ue[is]); 01769 } 01770 01771 for (unsigned int js = 0; js != n_side_dofs; ++js) 01772 { 01773 const unsigned int j = my_side_dofs[js]; 01774 const dof_id_type my_dof_g = my_dof_indices[j]; 01775 libmesh_assert_not_equal_to (my_dof_g, DofObject::invalid_id); 01776 01777 // Hunt for "constraining against myself" cases before 01778 // we bother creating a constraint row 01779 bool self_constraint = false; 01780 for (unsigned int is = 0; is != n_side_dofs; ++is) 01781 { 01782 const unsigned int i = neigh_side_dofs[is]; 01783 const dof_id_type their_dof_g = neigh_dof_indices[i]; 01784 libmesh_assert_not_equal_to (their_dof_g, DofObject::invalid_id); 01785 01786 if (their_dof_g == my_dof_g) 01787 { 01788 #ifndef NDEBUG 01789 const Real their_dof_value = Ue[is](js); 01790 libmesh_assert_less (std::abs(their_dof_value-1.), 01791 10*TOLERANCE); 01792 01793 for (unsigned int k = 0; k != n_side_dofs; ++k) 01794 libmesh_assert(k == is || 01795 std::abs(Ue[k](js)) < 01796 10*TOLERANCE); 01797 #endif 01798 01799 self_constraint = true; 01800 break; 01801 } 01802 } 01803 01804 if (self_constraint) 01805 continue; 01806 01807 DofConstraintRow* constraint_row; 01808 01809 // we may be running constraint methods concurretly 01810 // on multiple threads, so we need a lock to 01811 // ensure that this constraint is "ours" 01812 { 01813 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 01814 01815 if (dof_map.is_constrained_dof(my_dof_g)) 01816 continue; 01817 01818 constraint_row = &(constraints[my_dof_g].first); 01819 libmesh_assert(constraint_row->empty()); 01820 constraints[my_dof_g].second = 0; 01821 } 01822 01823 for (unsigned int is = 0; is != n_side_dofs; ++is) 01824 { 01825 const unsigned int i = neigh_side_dofs[is]; 01826 const dof_id_type their_dof_g = neigh_dof_indices[i]; 01827 libmesh_assert_not_equal_to (their_dof_g, DofObject::invalid_id); 01828 libmesh_assert_not_equal_to (their_dof_g, my_dof_g); 01829 01830 const Real their_dof_value = Ue[is](js); 01831 01832 if (std::abs(their_dof_value) < 10*TOLERANCE) 01833 continue; 01834 01835 constraint_row->insert(std::make_pair(their_dof_g, 01836 their_dof_value)); 01837 } 01838 } 01839 } 01840 // p refinement constraints: 01841 // constrain dofs shared between 01842 // active elements and neighbors with 01843 // lower polynomial degrees 01844 const unsigned int min_p_level = 01845 neigh->min_p_level_by_neighbor(elem, elem->p_level()); 01846 if (min_p_level < elem->p_level()) 01847 { 01848 // Adaptive p refinement of non-hierarchic bases will 01849 // require more coding 01850 libmesh_assert(my_fe->is_hierarchic()); 01851 dof_map.constrain_p_dofs(variable_number, elem, 01852 s, min_p_level); 01853 } 01854 } 01855 } 01856 01857 #endif // #ifdef LIBMESH_ENABLE_AMR 01858 01859 01860 01861 #ifdef LIBMESH_ENABLE_PERIODIC 01862 template <typename OutputType> 01863 void 01864 FEGenericBase<OutputType>:: 01865 compute_periodic_constraints (DofConstraints &constraints, 01866 DofMap &dof_map, 01867 const PeriodicBoundaries &boundaries, 01868 const MeshBase &mesh, 01869 const PointLocatorBase *point_locator, 01870 const unsigned int variable_number, 01871 const Elem* elem) 01872 { 01873 // Only bother if we truly have periodic boundaries 01874 if (boundaries.empty()) 01875 return; 01876 01877 libmesh_assert(elem); 01878 01879 // Only constrain active elements with this method 01880 if (!elem->active()) 01881 return; 01882 01883 const unsigned int Dim = elem->dim(); 01884 01885 // We need sys_number and variable_number for DofObject methods 01886 // later 01887 const unsigned int sys_number = dof_map.sys_number(); 01888 01889 const FEType& base_fe_type = dof_map.variable_type(variable_number); 01890 01891 // Construct FE objects for this element and its pseudo-neighbors. 01892 AutoPtr<FEGenericBase<OutputShape> > my_fe 01893 (FEGenericBase<OutputShape>::build(Dim, base_fe_type)); 01894 const FEContinuity cont = my_fe->get_continuity(); 01895 01896 // We don't need to constrain discontinuous elements 01897 if (cont == DISCONTINUOUS) 01898 return; 01899 libmesh_assert (cont == C_ZERO || cont == C_ONE); 01900 01901 // We'll use element size to generate relative tolerances later 01902 const Real primary_hmin = elem->hmin(); 01903 01904 AutoPtr<FEGenericBase<OutputShape> > neigh_fe 01905 (FEGenericBase<OutputShape>::build(Dim, base_fe_type)); 01906 01907 QGauss my_qface(Dim-1, base_fe_type.default_quadrature_order()); 01908 my_fe->attach_quadrature_rule (&my_qface); 01909 std::vector<Point> neigh_qface; 01910 01911 const std::vector<Real>& JxW = my_fe->get_JxW(); 01912 const std::vector<Point>& q_point = my_fe->get_xyz(); 01913 const std::vector<std::vector<OutputShape> >& phi = my_fe->get_phi(); 01914 const std::vector<std::vector<OutputShape> >& neigh_phi = 01915 neigh_fe->get_phi(); 01916 const std::vector<Point> *face_normals = NULL; 01917 const std::vector<std::vector<OutputGradient> > *dphi = NULL; 01918 const std::vector<std::vector<OutputGradient> > *neigh_dphi = NULL; 01919 std::vector<dof_id_type> my_dof_indices, neigh_dof_indices; 01920 std::vector<unsigned int> my_side_dofs, neigh_side_dofs; 01921 01922 if (cont != C_ZERO) 01923 { 01924 const std::vector<Point>& ref_face_normals = 01925 my_fe->get_normals(); 01926 face_normals = &ref_face_normals; 01927 const std::vector<std::vector<OutputGradient> >& ref_dphi = 01928 my_fe->get_dphi(); 01929 dphi = &ref_dphi; 01930 const std::vector<std::vector<OutputGradient> >& ref_neigh_dphi = 01931 neigh_fe->get_dphi(); 01932 neigh_dphi = &ref_neigh_dphi; 01933 } 01934 01935 DenseMatrix<Real> Ke; 01936 DenseVector<Real> Fe; 01937 std::vector<DenseVector<Real> > Ue; 01938 01939 // Look at the element faces. Check to see if we need to 01940 // build constraints. 01941 for (unsigned int s=0; s<elem->n_sides(); s++) 01942 { 01943 if (elem->neighbor(s)) 01944 continue; 01945 01946 const std::vector<boundary_id_type>& bc_ids = mesh.boundary_info->boundary_ids (elem, s); 01947 for (std::vector<boundary_id_type>::const_iterator id_it=bc_ids.begin(); id_it!=bc_ids.end(); ++id_it) 01948 { 01949 const boundary_id_type boundary_id = *id_it; 01950 const PeriodicBoundaryBase *periodic = boundaries.boundary(boundary_id); 01951 if (periodic && periodic->is_my_variable(variable_number)) 01952 { 01953 libmesh_assert(point_locator); 01954 01955 // Get pointers to the element's neighbor. 01956 const Elem* neigh = boundaries.neighbor(boundary_id, *point_locator, elem, s); 01957 01958 if (neigh == NULL) 01959 { 01960 libMesh::err << "PeriodicBoundaries point locator object returned NULL!" << std::endl; 01961 libmesh_error(); 01962 } 01963 01964 // periodic (and possibly h refinement) constraints: 01965 // constrain dofs shared between 01966 // this element and ones as coarse 01967 // as or coarser than this element. 01968 if (neigh->level() <= elem->level()) 01969 { 01970 unsigned int s_neigh = 01971 mesh.boundary_info->side_with_boundary_id (neigh, periodic->pairedboundary); 01972 libmesh_assert_not_equal_to (s_neigh, libMesh::invalid_uint); 01973 01974 #ifdef LIBMESH_ENABLE_AMR 01975 // Find the minimum p level; we build the h constraint 01976 // matrix with this and then constrain away all higher p 01977 // DoFs. 01978 libmesh_assert(neigh->active()); 01979 const unsigned int min_p_level = 01980 std::min(elem->p_level(), neigh->p_level()); 01981 01982 // we may need to make the FE objects reinit with the 01983 // minimum shared p_level 01984 // FIXME - I hate using const_cast<> and avoiding 01985 // accessor functions; there's got to be a 01986 // better way to do this! 01987 const unsigned int old_elem_level = elem->p_level(); 01988 if (old_elem_level != min_p_level) 01989 (const_cast<Elem *>(elem))->hack_p_level(min_p_level); 01990 const unsigned int old_neigh_level = neigh->p_level(); 01991 if (old_neigh_level != min_p_level) 01992 (const_cast<Elem *>(neigh))->hack_p_level(min_p_level); 01993 #endif // #ifdef LIBMESH_ENABLE_AMR 01994 01995 // We can do a projection with a single integration, 01996 // due to the assumption of nested finite element 01997 // subspaces. 01998 // FIXME: it might be more efficient to do nodes, 01999 // then edges, then side, to reduce the size of the 02000 // Cholesky factorization(s) 02001 my_fe->reinit(elem, s); 02002 02003 dof_map.dof_indices (elem, my_dof_indices, 02004 variable_number); 02005 dof_map.dof_indices (neigh, neigh_dof_indices, 02006 variable_number); 02007 02008 const unsigned int n_qp = my_qface.n_points(); 02009 02010 // Translate the quadrature points over to the 02011 // neighbor's boundary 02012 std::vector<Point> neigh_point(q_point.size()); 02013 for (unsigned int i=0; i != neigh_point.size(); ++i) 02014 neigh_point[i] = periodic->get_corresponding_pos(q_point[i]); 02015 02016 FEInterface::inverse_map (Dim, base_fe_type, neigh, 02017 neigh_point, neigh_qface); 02018 02019 neigh_fe->reinit(neigh, &neigh_qface); 02020 02021 // We're only concerned with DOFs whose values (and/or first 02022 // derivatives for C1 elements) are supported on side nodes 02023 FEInterface::dofs_on_side(elem, Dim, base_fe_type, s, my_side_dofs); 02024 FEInterface::dofs_on_side(neigh, Dim, base_fe_type, s_neigh, neigh_side_dofs); 02025 02026 // We're done with functions that examine Elem::p_level(), 02027 // so let's unhack those levels 02028 #ifdef LIBMESH_ENABLE_AMR 02029 if (elem->p_level() != old_elem_level) 02030 (const_cast<Elem *>(elem))->hack_p_level(old_elem_level); 02031 if (neigh->p_level() != old_neigh_level) 02032 (const_cast<Elem *>(neigh))->hack_p_level(old_neigh_level); 02033 #endif // #ifdef LIBMESH_ENABLE_AMR 02034 02035 const unsigned int n_side_dofs = 02036 libmesh_cast_int<unsigned int> 02037 (my_side_dofs.size()); 02038 libmesh_assert_equal_to (n_side_dofs, neigh_side_dofs.size()); 02039 02040 Ke.resize (n_side_dofs, n_side_dofs); 02041 Ue.resize(n_side_dofs); 02042 02043 // Form the projection matrix, (inner product of fine basis 02044 // functions against fine test functions) 02045 for (unsigned int is = 0; is != n_side_dofs; ++is) 02046 { 02047 const unsigned int i = my_side_dofs[is]; 02048 for (unsigned int js = 0; js != n_side_dofs; ++js) 02049 { 02050 const unsigned int j = my_side_dofs[js]; 02051 for (unsigned int qp = 0; qp != n_qp; ++qp) 02052 { 02053 Ke(is,js) += JxW[qp] * 02054 TensorTools::inner_product(phi[i][qp], 02055 phi[j][qp]); 02056 if (cont != C_ZERO) 02057 Ke(is,js) += JxW[qp] * 02058 TensorTools::inner_product((*dphi)[i][qp] * 02059 (*face_normals)[qp], 02060 (*dphi)[j][qp] * 02061 (*face_normals)[qp]); 02062 } 02063 } 02064 } 02065 02066 // Form the right hand sides, (inner product of coarse basis 02067 // functions against fine test functions) 02068 for (unsigned int is = 0; is != n_side_dofs; ++is) 02069 { 02070 const unsigned int i = neigh_side_dofs[is]; 02071 Fe.resize (n_side_dofs); 02072 for (unsigned int js = 0; js != n_side_dofs; ++js) 02073 { 02074 const unsigned int j = my_side_dofs[js]; 02075 for (unsigned int qp = 0; qp != n_qp; ++qp) 02076 { 02077 Fe(js) += JxW[qp] * 02078 TensorTools::inner_product(neigh_phi[i][qp], 02079 phi[j][qp]); 02080 if (cont != C_ZERO) 02081 Fe(js) += JxW[qp] * 02082 TensorTools::inner_product((*neigh_dphi)[i][qp] * 02083 (*face_normals)[qp], 02084 (*dphi)[j][qp] * 02085 (*face_normals)[qp]); 02086 } 02087 } 02088 Ke.cholesky_solve(Fe, Ue[is]); 02089 } 02090 02091 // Make sure we're not adding recursive constraints 02092 // due to the redundancy in the way we add periodic 02093 // boundary constraints 02094 // 02095 // In order for this to work while threaded or on 02096 // distributed meshes, we need a rigorous way to 02097 // avoid recursive constraints. Here it is: 02098 // 02099 // For vertex DoFs, if there is a "prior" element 02100 // (i.e. a coarser element or an equally refined 02101 // element with a lower id) on this boundary which 02102 // contains the vertex point, then we will avoid 02103 // generating constraints; the prior element (or 02104 // something prior to it) may do so. If we are the 02105 // most prior (or "primary") element on this 02106 // boundary sharing this point, then we look at the 02107 // boundary periodic to us, we find the primary 02108 // element there, and if that primary is coarser or 02109 // equal-but-lower-id, then our vertex dofs are 02110 // constrained in terms of that element. 02111 // 02112 // For edge DoFs, if there is a coarser element 02113 // on this boundary sharing this edge, then we will 02114 // avoid generating constraints (we will be 02115 // constrained indirectly via AMR constraints 02116 // connecting us to the coarser element's DoFs). If 02117 // we are the coarsest element sharing this edge, 02118 // then we generate constraints if and only if we 02119 // are finer than the coarsest element on the 02120 // boundary periodic to us sharing the corresponding 02121 // periodic edge, or if we are at equal level but 02122 // our edge nodes have higher ids than the periodic 02123 // edge nodes (sorted from highest to lowest, then 02124 // compared lexicographically) 02125 // 02126 // For face DoFs, we generate constraints if we are 02127 // finer than our periodic neighbor, or if we are at 02128 // equal level but our element id is higher than its 02129 // element id. 02130 // 02131 // If the primary neighbor is also the current elem 02132 // (a 1-element-thick mesh) then we choose which 02133 // vertex dofs to constrain via lexicographic 02134 // ordering on point locations 02135 02136 // FIXME: This code doesn't yet properly handle 02137 // cases where multiple different periodic BCs 02138 // intersect. 02139 std::set<dof_id_type> my_constrained_dofs; 02140 02141 for (unsigned int n = 0; n != elem->n_nodes(); ++n) 02142 { 02143 if (!elem->is_node_on_side(n,s)) 02144 continue; 02145 02146 const Node* my_node = elem->get_node(n); 02147 02148 if (elem->is_vertex(n)) 02149 { 02150 // Find all boundary ids that include this 02151 // point and have periodic boundary 02152 // conditions for this variable 02153 std::set<boundary_id_type> point_bcids; 02154 02155 for (unsigned int new_s = 0; new_s != 02156 elem->n_sides(); ++new_s) 02157 { 02158 if (!elem->is_node_on_side(n,new_s)) 02159 continue; 02160 02161 const std::vector<boundary_id_type> 02162 new_bc_ids = mesh.boundary_info->boundary_ids (elem, s); 02163 for (std::vector<boundary_id_type>::const_iterator 02164 new_id_it=new_bc_ids.begin(); new_id_it!=new_bc_ids.end(); ++new_id_it) 02165 { 02166 const boundary_id_type new_boundary_id = *new_id_it; 02167 const PeriodicBoundaryBase *new_periodic = boundaries.boundary(new_boundary_id); 02168 if (new_periodic && new_periodic->is_my_variable(variable_number)) 02169 { 02170 point_bcids.insert(new_boundary_id); 02171 } 02172 } 02173 } 02174 02175 // See if this vertex has point neighbors to 02176 // defer to 02177 if (primary_boundary_point_neighbor 02178 (elem, *my_node, *mesh.boundary_info, point_bcids) != elem) 02179 continue; 02180 02181 // Find the complementary boundary id set 02182 std::set<boundary_id_type> point_pairedids; 02183 for (std::set<boundary_id_type>::const_iterator i = 02184 point_bcids.begin(); i != point_bcids.end(); ++i) 02185 { 02186 const boundary_id_type new_boundary_id = *i; 02187 const PeriodicBoundaryBase *new_periodic = boundaries.boundary(new_boundary_id); 02188 point_pairedids.insert(new_periodic->pairedboundary); 02189 } 02190 02191 // What do we want to constrain against? 02192 const Elem* primary_elem = NULL; 02193 const Elem* main_neigh = NULL; 02194 Point main_pt = *my_node, 02195 primary_pt = *my_node; 02196 02197 for (std::set<boundary_id_type>::const_iterator i = 02198 point_bcids.begin(); i != point_bcids.end(); ++i) 02199 { 02200 // Find the corresponding periodic point and 02201 // its primary neighbor 02202 const boundary_id_type new_boundary_id = *i; 02203 const PeriodicBoundaryBase *new_periodic = boundaries.boundary(new_boundary_id); 02204 02205 const Point neigh_pt = 02206 new_periodic->get_corresponding_pos(*my_node); 02207 02208 // If the point is getting constrained 02209 // to itself by this PBC then we don't 02210 // generate any constraints 02211 if (neigh_pt.absolute_fuzzy_equals 02212 (*my_node, primary_hmin*TOLERANCE)) 02213 continue; 02214 02215 // Otherwise we'll have a constraint in 02216 // one direction or another 02217 if (!primary_elem) 02218 primary_elem = elem; 02219 02220 const Elem *primary_neigh = primary_boundary_point_neighbor 02221 (neigh, neigh_pt, *mesh.boundary_info, 02222 point_pairedids); 02223 02224 libmesh_assert(primary_neigh); 02225 02226 if (new_boundary_id == boundary_id) 02227 { 02228 main_neigh = primary_neigh; 02229 main_pt = neigh_pt; 02230 } 02231 02232 // Finer elements will get constrained in 02233 // terms of coarser neighbors, not the 02234 // other way around 02235 if ((primary_neigh->level() > primary_elem->level()) || 02236 02237 // For equal-level elements, the one with 02238 // higher id gets constrained in terms of 02239 // the one with lower id 02240 (primary_neigh->level() == primary_elem->level() && 02241 primary_neigh->id() > primary_elem->id()) || 02242 02243 // On a one-element-thick mesh, we compare 02244 // points to see what side gets constrained 02245 (primary_neigh == primary_elem && 02246 (neigh_pt > primary_pt))) 02247 continue; 02248 02249 primary_elem = primary_neigh; 02250 primary_pt = neigh_pt; 02251 } 02252 02253 if (!primary_elem || 02254 primary_elem != main_neigh || 02255 primary_pt != main_pt) 02256 continue; 02257 } 02258 else if (elem->is_edge(n)) 02259 { 02260 // Find which edge we're on 02261 unsigned int e=0; 02262 for (; e != elem->n_edges(); ++e) 02263 { 02264 if (elem->is_node_on_edge(n,e)) 02265 break; 02266 } 02267 libmesh_assert_less (e, elem->n_edges()); 02268 02269 // Find the edge end nodes 02270 Node *e1 = NULL, 02271 *e2 = NULL; 02272 for (unsigned int nn = 0; nn != elem->n_nodes(); ++nn) 02273 { 02274 if (nn == n) 02275 continue; 02276 02277 if (elem->is_node_on_edge(nn, e)) 02278 { 02279 if (e1 == NULL) 02280 { 02281 e1 = elem->get_node(nn); 02282 } 02283 else 02284 { 02285 e2 = elem->get_node(nn); 02286 break; 02287 } 02288 } 02289 } 02290 libmesh_assert (e1 && e2); 02291 02292 // Find all boundary ids that include this 02293 // edge and have periodic boundary 02294 // conditions for this variable 02295 std::set<boundary_id_type> edge_bcids; 02296 02297 for (unsigned int new_s = 0; new_s != 02298 elem->n_sides(); ++new_s) 02299 { 02300 if (!elem->is_node_on_side(n,new_s)) 02301 continue; 02302 02303 const std::vector<boundary_id_type>& 02304 new_bc_ids = mesh.boundary_info->boundary_ids (elem, s); 02305 for (std::vector<boundary_id_type>::const_iterator 02306 new_id_it=new_bc_ids.begin(); new_id_it!=new_bc_ids.end(); ++new_id_it) 02307 { 02308 const boundary_id_type new_boundary_id = *new_id_it; 02309 const PeriodicBoundaryBase *new_periodic = boundaries.boundary(new_boundary_id); 02310 if (new_periodic && new_periodic->is_my_variable(variable_number)) 02311 { 02312 edge_bcids.insert(new_boundary_id); 02313 } 02314 } 02315 } 02316 02317 02318 // See if this edge has neighbors to defer to 02319 if (primary_boundary_edge_neighbor 02320 (elem, *e1, *e2, *mesh.boundary_info, edge_bcids) != elem) 02321 continue; 02322 02323 // Find the complementary boundary id set 02324 std::set<boundary_id_type> edge_pairedids; 02325 for (std::set<boundary_id_type>::const_iterator i = 02326 edge_bcids.begin(); i != edge_bcids.end(); ++i) 02327 { 02328 const boundary_id_type new_boundary_id = *i; 02329 const PeriodicBoundaryBase *new_periodic = boundaries.boundary(new_boundary_id); 02330 edge_pairedids.insert(new_periodic->pairedboundary); 02331 } 02332 02333 // What do we want to constrain against? 02334 const Elem* primary_elem = NULL; 02335 const Elem* main_neigh = NULL; 02336 Point main_pt1 = *e1, 02337 main_pt2 = *e2, 02338 primary_pt1 = *e1, 02339 primary_pt2 = *e2; 02340 02341 for (std::set<boundary_id_type>::const_iterator i = 02342 edge_bcids.begin(); i != edge_bcids.end(); ++i) 02343 { 02344 // Find the corresponding periodic edge and 02345 // its primary neighbor 02346 const boundary_id_type new_boundary_id = *i; 02347 const PeriodicBoundaryBase *new_periodic = boundaries.boundary(new_boundary_id); 02348 02349 Point neigh_pt1 = new_periodic->get_corresponding_pos(*e1), 02350 neigh_pt2 = new_periodic->get_corresponding_pos(*e2); 02351 02352 // If the edge is getting constrained 02353 // to itself by this PBC then we don't 02354 // generate any constraints 02355 if (neigh_pt1.absolute_fuzzy_equals 02356 (*e1, primary_hmin*TOLERANCE) && 02357 neigh_pt2.absolute_fuzzy_equals 02358 (*e2, primary_hmin*TOLERANCE)) 02359 continue; 02360 02361 // Otherwise we'll have a constraint in 02362 // one direction or another 02363 if (!primary_elem) 02364 primary_elem = elem; 02365 02366 const Elem *primary_neigh = primary_boundary_edge_neighbor 02367 (neigh, neigh_pt1, neigh_pt2, *mesh.boundary_info, 02368 edge_pairedids); 02369 02370 libmesh_assert(primary_neigh); 02371 02372 if (new_boundary_id == boundary_id) 02373 { 02374 main_neigh = primary_neigh; 02375 main_pt1 = neigh_pt1; 02376 main_pt2 = neigh_pt2; 02377 } 02378 02379 // If we have a one-element thick mesh, 02380 // we'll need to sort our points to get a 02381 // consistent ordering rule 02382 // 02383 // Use >= in this test to make sure that, 02384 // for angular constraints, no node gets 02385 // constrained to itself. 02386 if (primary_neigh == primary_elem) 02387 { 02388 if (primary_pt1 > primary_pt2) 02389 std::swap(primary_pt1, primary_pt2); 02390 if (neigh_pt1 > neigh_pt2) 02391 std::swap(neigh_pt1, neigh_pt2); 02392 02393 if (neigh_pt2 >= primary_pt2) 02394 continue; 02395 } 02396 02397 // Otherwise: 02398 // Finer elements will get constrained in 02399 // terms of coarser ones, not the other way 02400 // around 02401 if ((primary_neigh->level() > primary_elem->level()) || 02402 02403 // For equal-level elements, the one with 02404 // higher id gets constrained in terms of 02405 // the one with lower id 02406 (primary_neigh->level() == primary_elem->level() && 02407 primary_neigh->id() > primary_elem->id())) 02408 continue; 02409 02410 primary_elem = primary_neigh; 02411 primary_pt1 = neigh_pt1; 02412 primary_pt2 = neigh_pt2; 02413 } 02414 02415 if (!primary_elem || 02416 primary_elem != main_neigh || 02417 primary_pt1 != main_pt1 || 02418 primary_pt2 != main_pt2) 02419 continue; 02420 } 02421 else if (elem->is_face(n)) 02422 { 02423 // If we have a one-element thick mesh, 02424 // use the ordering of the face node and its 02425 // periodic counterpart to determine what 02426 // gets constrained 02427 if (neigh == elem) 02428 { 02429 const Point neigh_pt = 02430 periodic->get_corresponding_pos(*my_node); 02431 if (neigh_pt > *my_node) 02432 continue; 02433 } 02434 02435 // Otherwise: 02436 // Finer elements will get constrained in 02437 // terms of coarser ones, not the other way 02438 // around 02439 if ((neigh->level() > elem->level()) || 02440 02441 // For equal-level elements, the one with 02442 // higher id gets constrained in terms of 02443 // the one with lower id 02444 (neigh->level() == elem->level() && 02445 neigh->id() > elem->id())) 02446 continue; 02447 } 02448 02449 // If we made it here without hitting a continue 02450 // statement, then we're at a node whose dofs 02451 // should be constrained by this element's 02452 // calculations. 02453 const unsigned int n_comp = 02454 my_node->n_comp(sys_number, variable_number); 02455 02456 for (unsigned int i=0; i != n_comp; ++i) 02457 my_constrained_dofs.insert 02458 (my_node->dof_number 02459 (sys_number, variable_number, i)); 02460 } 02461 02462 // FIXME: old code for disambiguating periodic BCs: 02463 // this is not threadsafe nor safe to run on a 02464 // non-serialized mesh. 02465 /* 02466 std::vector<bool> recursive_constraint(n_side_dofs, false); 02467 02468 for (unsigned int is = 0; is != n_side_dofs; ++is) 02469 { 02470 const unsigned int i = neigh_side_dofs[is]; 02471 const dof_id_type their_dof_g = neigh_dof_indices[i]; 02472 libmesh_assert_not_equal_to (their_dof_g, DofObject::invalid_id); 02473 02474 { 02475 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 02476 02477 if (!dof_map.is_constrained_dof(their_dof_g)) 02478 continue; 02479 } 02480 02481 DofConstraintRow& their_constraint_row = 02482 constraints[their_dof_g].first; 02483 02484 for (unsigned int js = 0; js != n_side_dofs; ++js) 02485 { 02486 const unsigned int j = my_side_dofs[js]; 02487 const dof_id_type my_dof_g = my_dof_indices[j]; 02488 libmesh_assert_not_equal_to (my_dof_g, DofObject::invalid_id); 02489 02490 if (their_constraint_row.count(my_dof_g)) 02491 recursive_constraint[js] = true; 02492 } 02493 } 02494 */ 02495 02496 for (unsigned int js = 0; js != n_side_dofs; ++js) 02497 { 02498 // FIXME: old code path 02499 // if (recursive_constraint[js]) 02500 // continue; 02501 02502 const unsigned int j = my_side_dofs[js]; 02503 const dof_id_type my_dof_g = my_dof_indices[j]; 02504 libmesh_assert_not_equal_to (my_dof_g, DofObject::invalid_id); 02505 02506 // FIXME: new code path 02507 if (!my_constrained_dofs.count(my_dof_g)) 02508 continue; 02509 02510 DofConstraintRow* constraint_row; 02511 02512 // we may be running constraint methods concurretly 02513 // on multiple threads, so we need a lock to 02514 // ensure that this constraint is "ours" 02515 { 02516 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 02517 02518 if (dof_map.is_constrained_dof(my_dof_g)) 02519 continue; 02520 02521 constraint_row = &(constraints[my_dof_g].first); 02522 libmesh_assert(constraint_row->empty()); 02523 constraints[my_dof_g].second = 0; 02524 } 02525 02526 for (unsigned int is = 0; is != n_side_dofs; ++is) 02527 { 02528 const unsigned int i = neigh_side_dofs[is]; 02529 const dof_id_type their_dof_g = neigh_dof_indices[i]; 02530 libmesh_assert_not_equal_to (their_dof_g, DofObject::invalid_id); 02531 02532 // Periodic constraints should never be 02533 // self-constraints 02534 // libmesh_assert_not_equal_to (their_dof_g, my_dof_g); 02535 02536 const Real their_dof_value = Ue[is](js); 02537 02538 if (their_dof_g == my_dof_g) 02539 { 02540 libmesh_assert_less (std::abs(their_dof_value-1.), 1.e-5); 02541 for (unsigned int k = 0; k != n_side_dofs; ++k) 02542 libmesh_assert(k == is || std::abs(Ue[k](js)) < 1.e-5); 02543 continue; 02544 } 02545 02546 if (std::abs(their_dof_value) < 10*TOLERANCE) 02547 continue; 02548 02549 constraint_row->insert(std::make_pair(their_dof_g, 02550 their_dof_value)); 02551 } 02552 } 02553 } 02554 // p refinement constraints: 02555 // constrain dofs shared between 02556 // active elements and neighbors with 02557 // lower polynomial degrees 02558 #ifdef LIBMESH_ENABLE_AMR 02559 const unsigned int min_p_level = 02560 neigh->min_p_level_by_neighbor(elem, elem->p_level()); 02561 if (min_p_level < elem->p_level()) 02562 { 02563 // Adaptive p refinement of non-hierarchic bases will 02564 // require more coding 02565 libmesh_assert(my_fe->is_hierarchic()); 02566 dof_map.constrain_p_dofs(variable_number, elem, 02567 s, min_p_level); 02568 } 02569 #endif // #ifdef LIBMESH_ENABLE_AMR 02570 } 02571 } 02572 } 02573 } 02574 02575 #endif // LIBMESH_ENABLE_PERIODIC 02576 02577 // ------------------------------------------------------------ 02578 // Explicit instantiations 02579 template class FEGenericBase<Real>; 02580 template class FEGenericBase<RealGradient>; 02581 02582 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: