inf_fe_map.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/libmesh_config.h" 00022 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00023 #include "libmesh/inf_fe.h" 00024 #include "libmesh/fe.h" 00025 #include "libmesh/elem.h" 00026 #include "libmesh/inf_fe_macro.h" 00027 #include "libmesh/libmesh_logging.h" 00028 00029 namespace libMesh 00030 { 00031 00032 00033 00034 // ------------------------------------------------------------ 00035 // InfFE static class members concerned with coordinate 00036 // mapping 00037 00038 00039 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00040 Point InfFE<Dim,T_radial,T_map>::map (const Elem* inf_elem, 00041 const Point& reference_point) 00042 { 00043 libmesh_assert(inf_elem); 00044 libmesh_assert_not_equal_to (Dim, 0); 00045 00046 AutoPtr<Elem> base_elem (Base::build_elem (inf_elem)); 00047 00048 const Order radial_mapping_order (Radial::mapping_order()); 00049 const Real v (reference_point(Dim-1)); 00050 00051 // map in the base face 00052 Point base_point; 00053 switch (Dim) 00054 { 00055 case 1: 00056 base_point = inf_elem->point(0); 00057 break; 00058 case 2: 00059 base_point = FE<1,LAGRANGE>::map (base_elem.get(), reference_point); 00060 break; 00061 case 3: 00062 base_point = FE<2,LAGRANGE>::map (base_elem.get(), reference_point); 00063 break; 00064 #ifdef DEBUG 00065 default: 00066 libmesh_error(); 00067 #endif 00068 } 00069 00070 00071 // map in the outer node face not necessary. Simply 00072 // compute the outer_point = base_point + (base_point-origin) 00073 const Point outer_point (base_point * 2. - inf_elem->origin()); 00074 00075 Point p; 00076 00077 // there are only two mapping shapes in radial direction 00078 p.add_scaled (base_point, InfFE<Dim,INFINITE_MAP,T_map>::eval (v, radial_mapping_order, 0)); 00079 p.add_scaled (outer_point, InfFE<Dim,INFINITE_MAP,T_map>::eval (v, radial_mapping_order, 1)); 00080 00081 return p; 00082 } 00083 00084 00085 00086 00087 00088 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00089 Point InfFE<Dim,T_radial,T_map>::inverse_map (const Elem* inf_elem, 00090 const Point& physical_point, 00091 const Real tolerance, 00092 const bool secure, 00093 const bool interpolated) 00094 { 00095 libmesh_assert(inf_elem); 00096 libmesh_assert_greater_equal (tolerance, 0.); 00097 00098 00099 // Start logging the map inversion. 00100 START_LOG("inverse_map()", "InfFE"); 00101 00102 00103 // 1.) 00104 // build a base element to do the map inversion in the base face 00105 AutoPtr<Elem> base_elem (Base::build_elem (inf_elem)); 00106 00107 const Order base_mapping_order (base_elem->default_order()); 00108 const ElemType base_mapping_elem_type (base_elem->type()); 00109 const unsigned int n_base_mapping_sf = Base::n_base_mapping_sf (base_mapping_elem_type, 00110 base_mapping_order); 00111 00112 const ElemType inf_elem_type = inf_elem->type(); 00113 if (inf_elem_type != INFHEX8 && 00114 inf_elem_type != INFPRISM6) 00115 { 00116 libMesh::err << "ERROR: InfFE::inverse_map is currently implemented only for\n" 00117 << " infinite elments of type InfHex8 and InfPrism6." << std::endl; 00118 libmesh_error(); 00119 } 00120 00121 00122 // 2.) 00123 // just like in FE<Dim-1,LAGRANGE>::inverse_map(): compute 00124 // the local coordinates, but only in the base element. 00125 // The radial part can then be computed directly later on. 00126 00127 00128 // How much did the point on the reference 00129 // element change by in this Newton step? 00130 Real inverse_map_error = 0.; 00131 00132 00133 // The point on the reference element. This is 00134 // the "initial guess" for Newton's method. The 00135 // centroid seems like a good idea, but computing 00136 // it is a little more intensive than, say taking 00137 // the zero point. 00138 // 00139 // Convergence should be insensitive of this choice 00140 // for "good" elements. 00141 Point p; // the zero point. No computation required 00142 00143 00144 00145 // Now find the intersection of a plane represented by the base 00146 // element nodes and the line given by the origin of the infinite 00147 // element and the physical point. 00148 Point intersection; 00149 00150 // the origin of the infinite lement 00151 const Point o = inf_elem->origin(); 00152 00153 switch (Dim) 00154 { 00155 00156 // unnecessary for 1D 00157 case 1: 00158 { 00159 break; 00160 } 00161 00162 case 2: 00163 { 00164 libMesh::err << "ERROR: InfFE::inverse_map is not yet implemented" 00165 << " in 2d" << std::endl; 00166 libmesh_error(); 00167 break; 00168 } 00169 00170 00171 case 3: 00172 { 00173 // references to the nodal points of the base element 00174 const Point& p0 = base_elem->point(0); 00175 const Point& p1 = base_elem->point(1); 00176 const Point& p2 = base_elem->point(2); 00177 00178 // a reference to the physical point 00179 const Point& fp = physical_point; 00180 00181 // The intersection of the plane and the line is given by 00182 // can be computed solving a linear 3x3 system 00183 // a*({p1}-{p0})+b*({p2}-{p0})-c*({fp}-{o})={fp}-{p0}. 00184 const Real c_factor = -(p1(0)*fp(1)*p0(2)-p1(0)*fp(2)*p0(1) 00185 +fp(0)*p1(2)*p0(1)-p0(0)*fp(1)*p1(2) 00186 +p0(0)*fp(2)*p1(1)+p2(0)*fp(2)*p0(1) 00187 -p2(0)*fp(1)*p0(2)-fp(0)*p2(2)*p0(1) 00188 +fp(0)*p0(2)*p2(1)+p0(0)*fp(1)*p2(2) 00189 -p0(0)*fp(2)*p2(1)-fp(0)*p0(2)*p1(1) 00190 +p0(2)*p2(0)*p1(1)-p0(1)*p2(0)*p1(2) 00191 -fp(0)*p1(2)*p2(1)+p2(1)*p0(0)*p1(2) 00192 -p2(0)*fp(2)*p1(1)-p1(0)*fp(1)*p2(2) 00193 +p2(2)*p1(0)*p0(1)+p1(0)*fp(2)*p2(1) 00194 -p0(2)*p1(0)*p2(1)-p2(2)*p0(0)*p1(1) 00195 +fp(0)*p2(2)*p1(1)+p2(0)*fp(1)*p1(2))/ 00196 (fp(0)*p1(2)*p0(1)-p1(0)*fp(2)*p0(1) 00197 +p1(0)*fp(1)*p0(2)-p1(0)*o(1)*p0(2) 00198 +o(0)*p2(2)*p0(1)-p0(0)*fp(2)*p2(1) 00199 +p1(0)*o(1)*p2(2)+fp(0)*p0(2)*p2(1) 00200 -fp(0)*p1(2)*p2(1)-p0(0)*o(1)*p2(2) 00201 +p0(0)*fp(1)*p2(2)-o(0)*p0(2)*p2(1) 00202 +o(0)*p1(2)*p2(1)-p2(0)*fp(2)*p1(1) 00203 +fp(0)*p2(2)*p1(1)-p2(0)*fp(1)*p0(2) 00204 -o(2)*p0(0)*p1(1)-fp(0)*p0(2)*p1(1) 00205 +p0(0)*o(1)*p1(2)+p0(0)*fp(2)*p1(1) 00206 -p0(0)*fp(1)*p1(2)-o(0)*p1(2)*p0(1) 00207 -p2(0)*o(1)*p1(2)-o(2)*p2(0)*p0(1) 00208 -o(2)*p1(0)*p2(1)+o(2)*p0(0)*p2(1) 00209 -fp(0)*p2(2)*p0(1)+o(2)*p2(0)*p1(1) 00210 +p2(0)*o(1)*p0(2)+p2(0)*fp(1)*p1(2) 00211 +p2(0)*fp(2)*p0(1)-p1(0)*fp(1)*p2(2) 00212 +p1(0)*fp(2)*p2(1)-o(0)*p2(2)*p1(1) 00213 +o(2)*p1(0)*p0(1)+o(0)*p0(2)*p1(1)); 00214 00215 00216 // Compute the intersection with 00217 // {intersection} = {fp} + c*({fp}-{o}). 00218 intersection.add_scaled(fp,1.+c_factor); 00219 intersection.add_scaled(o,-c_factor); 00220 00221 break; 00222 } 00223 } 00224 00228 unsigned int cnt = 0; 00229 00230 00234 do 00235 { 00236 00237 // Increment in current iterate \p p, will be computed. 00238 // Automatically initialized to all zero. Note that 00239 // in 3D, actually only the first two entries are 00240 // filled by the inverse map, and in 2D only the first. 00241 Point dp; 00242 00243 // The form of the map and how we invert it depends 00244 // on the dimension that we are in. 00245 switch (Dim) 00246 { 00247 00248 //------------------------------------------------------------------ 00249 // 1D infinite element - no map inversion necessary 00250 case 1: 00251 { 00252 break; 00253 } 00254 00255 00256 //------------------------------------------------------------------ 00257 // 2D infinite element - 1D map inversion 00258 // 00259 // In this iteration scheme only search for the local coordinate 00260 // in xi direction. Once xi is determined, the radial coordinate eta is 00261 // uniquely determined, and there is no need to iterate in that direction. 00262 case 2: 00263 { 00264 00265 // Where our current iterate \p p maps to. 00266 const Point physical_guess = FE<1,LAGRANGE>::map (base_elem.get(), p); 00267 00268 00269 // How far our current iterate is from the actual point. 00270 const Point delta = physical_point - physical_guess; 00271 00272 00273 const Point dxi = FE<1,LAGRANGE>::map_xi (base_elem.get(), p); 00274 00275 00276 // For details on Newton's method see fe_map.C 00277 const Real G = dxi*dxi; 00278 00279 if (secure) 00280 libmesh_assert_greater (G, 0.); 00281 00282 const Real Ginv = 1./G; 00283 00284 const Real dxidelta = dxi*delta; 00285 00286 // compute only the first coordinate 00287 dp(0) = Ginv*dxidelta; 00288 00289 break; 00290 } 00291 00292 00293 00294 //------------------------------------------------------------------ 00295 // 3D infinite element - 2D map inversion 00296 // 00297 // In this iteration scheme only search for the local coordinates 00298 // in xi and eta direction. Once xi, eta are determined, the radial 00299 // coordinate zeta may directly computed. 00300 case 3: 00301 { 00302 00303 00304 // Where our current iterate \p p maps to. 00305 const Point physical_guess = FE<2,LAGRANGE>::map (base_elem.get(), p); 00306 00307 // How far our current iterate is from the actual point. 00308 // const Point delta = physical_point - physical_guess; 00309 const Point delta = intersection - physical_guess; 00310 00311 const Point dxi = FE<2,LAGRANGE>::map_xi (base_elem.get(), p); 00312 const Point deta = FE<2,LAGRANGE>::map_eta (base_elem.get(), p); 00313 00314 00315 // For details on Newton's method see fe_map.C 00316 const Real 00317 G11 = dxi*dxi, G12 = dxi*deta, 00318 G21 = dxi*deta, G22 = deta*deta; 00319 00320 00321 const Real det = (G11*G22 - G12*G21); 00322 00323 if (secure) 00324 { 00325 libmesh_assert_greater (det, 0.); 00326 libmesh_assert_greater (std::abs(det), 1.e-10); 00327 } 00328 00329 const Real inv_det = 1./det; 00330 00331 const Real 00332 Ginv11 = G22*inv_det, 00333 Ginv12 = -G12*inv_det, 00334 00335 Ginv21 = -G21*inv_det, 00336 Ginv22 = G11*inv_det; 00337 00338 00339 const Real dxidelta = dxi*delta; 00340 const Real detadelta = deta*delta; 00341 00342 // compute only the first two coordinates. 00343 dp(0) = (Ginv11*dxidelta + Ginv12*detadelta); 00344 dp(1) = (Ginv21*dxidelta + Ginv22*detadelta); 00345 00346 break; 00347 } 00348 00349 00350 00351 // Some other dimension? 00352 default: 00353 libmesh_error(); 00354 } // end switch(Dim), dp now computed 00355 00356 00357 00358 // determine the error in computing the local coordinates 00359 // in the base: ||P_n+1 - P_n|| 00360 inverse_map_error = dp.size(); 00361 00362 00363 // P_n+1 = P_n + dp 00364 p.add (dp); 00365 00366 00367 // Increment the iteration count. 00368 cnt++; 00369 00370 00371 // Watch for divergence of Newton's 00372 // method. 00373 if (cnt > 10) 00374 { 00375 if (secure || !secure) 00376 { 00377 libmesh_here(); 00378 { 00379 libMesh::err << "WARNING: Newton scheme has not converged in " 00380 << cnt << " iterations:\n" 00381 << " physical_point=" 00382 << physical_point 00383 << " dp=" 00384 << dp 00385 << " p=" 00386 << p 00387 << " error=" << inverse_map_error 00388 << std::endl; 00389 } 00390 } 00391 00392 if (cnt > 20) 00393 { 00394 libMesh::err << "ERROR: Newton scheme FAILED to converge in " 00395 << cnt << " iterations!" << std::endl; 00396 libmesh_error(); 00397 } 00398 00399 // else 00400 // { 00401 // break; 00402 // } 00403 } 00404 } 00405 while (inverse_map_error > tolerance); 00406 00407 00408 00409 // 4. 00410 // 00411 // Now that we have the local coordinates in the base, 00412 // compute the interpolated radial distance a(s,t) \p a_interpolated 00413 if (interpolated) 00414 switch (Dim) 00415 { 00416 case 1: 00417 { 00418 Real a_interpolated = Point( inf_elem->point(0) 00419 - inf_elem->point(n_base_mapping_sf) ).size(); 00420 00421 p(0) = 1. - 2*a_interpolated/physical_point(0); 00422 00423 #ifdef DEBUG 00424 // the radial distance should always be >= -1. 00425 00426 if (p(0)+1 < tolerance) 00427 { 00428 libmesh_here(); 00429 libMesh::err << "WARNING: radial distance p(0) is " 00430 << p(0) 00431 << std::endl; 00432 } 00433 #endif 00434 00435 break; 00436 } 00437 00438 00439 case 2: 00440 { 00441 Real a_interpolated = 0.; 00442 00443 // the distance between the origin and the physical point 00444 const Real fp_o_dist = Point(o-physical_point).size(); 00445 00446 for (unsigned int i=0; i<n_base_mapping_sf; i++) 00447 { 00448 // the radial distance of the i-th base mapping point 00449 const Real dist_i = Point( inf_elem->point(i) 00450 - inf_elem->point(i+n_base_mapping_sf) ).size(); 00451 // weight with the corresponding shape function 00452 a_interpolated += dist_i * FE<1,LAGRANGE>::shape(base_mapping_elem_type, 00453 base_mapping_order, 00454 i, 00455 p); 00456 } 00457 00458 p(1) = 1. - 2*a_interpolated/fp_o_dist; 00459 00460 #ifdef DEBUG 00461 // the radial distance should always be >= -1. 00462 00463 // if (p(1)+1 < tolerance) 00464 // { 00465 // libmesh_here(); 00466 // libMesh::err << "WARNING: radial distance p(1) is " 00467 // << p(1) 00468 // << std::endl; 00469 // } 00470 #endif 00471 00472 break; 00473 } 00474 00475 00476 case 3: 00477 { 00478 Real a_interpolated = 0.; 00479 00480 00481 // the distance between the origin and the physical point 00482 const Real fp_o_dist = Point(o-physical_point).size(); 00483 00484 for (unsigned int i=0; i<n_base_mapping_sf; i++) 00485 { 00486 // the radial distance of the i-th base mapping point 00487 const Real dist_i = Point( inf_elem->point(i) 00488 - inf_elem->point(i+n_base_mapping_sf) ).size(); 00489 00490 // weight with the corresponding shape function 00491 a_interpolated += dist_i * FE<2,LAGRANGE>::shape(base_mapping_elem_type, 00492 base_mapping_order, 00493 i, 00494 p); 00495 00496 } 00497 00498 p(2) = 1. - 2*a_interpolated/fp_o_dist; 00499 00500 #ifdef DEBUG 00501 00502 00503 // the radial distance should always be >= -1. 00504 00505 // if (p(2)+1 < tolerance) 00506 // { 00507 // libmesh_here(); 00508 // libMesh::err << "WARNING: radial distance p(2) is " 00509 // << p(2) 00510 // << std::endl; 00511 // } 00512 #endif 00513 00514 break; 00515 } 00516 00517 default: 00518 libmesh_error(); 00519 } // end switch(Dim), p fully computed, including radial part 00520 00521 // if we do not want the interpolated distance, then 00522 // use newton iteration to get the actual distance 00523 else 00524 { 00525 // distance from the physical point to the ifem origin 00526 const Real fp_o_dist = Point(o-physical_point).size(); 00527 00528 // the distance from the intersection on the 00529 // base to the origin 00530 const Real a_dist = intersection.size(); 00531 00532 // element coordinate in radial direction 00533 // here our first guess is 0. 00534 Real v = 0.; 00535 00536 // the order of the radial mapping 00537 const Order radial_mapping_order (Radial::mapping_order()); 00538 00539 unsigned int cnt = 0; 00540 inverse_map_error = 0.; 00541 00542 // Newton iteration in 1-D 00543 do 00544 { 00545 // the mapping in radial direction 00546 // note that we only have two mapping functions in 00547 // radial direction 00548 const Real r = a_dist * InfFE<Dim,INFINITE_MAP,T_map>::eval (v, radial_mapping_order, 0) 00549 + 2. * a_dist * InfFE<Dim,INFINITE_MAP,T_map>::eval (v, radial_mapping_order, 1); 00550 00551 const Real dr = a_dist * InfFE<Dim,INFINITE_MAP,T_map>::eval_deriv (v, radial_mapping_order, 0) 00552 + 2. * a_dist * InfFE<Dim,INFINITE_MAP,T_map>::eval_deriv (v, radial_mapping_order, 1); 00553 00554 const Real G = dr*dr; 00555 const Real Ginv = 1./G; 00556 00557 const Real delta = fp_o_dist - r; 00558 const Real drdelta = dr*delta; 00559 00560 Real dp = Ginv*drdelta; 00561 00562 // update the radial coordinate 00563 v += dp; 00564 00565 // note that v should be smaller than 1, 00566 // since radial mapping function tends to infinity 00567 if (v >= 1.) 00568 v = .9999; 00569 00570 inverse_map_error = std::fabs(dp); 00571 00572 // increment iteration count 00573 cnt ++; 00574 if (cnt > 20) 00575 { 00576 00577 libMesh::err << "ERROR: 1D Newton scheme FAILED to converge" 00578 << std::endl; 00579 00580 libmesh_error(); 00581 } 00582 00583 00584 } 00585 while (inverse_map_error > tolerance); 00586 00587 switch (Dim) 00588 { 00589 case 1: 00590 { 00591 p(0) = v; 00592 break; 00593 } 00594 case 2: 00595 { 00596 p(1) = v; 00597 break; 00598 } 00599 case 3: 00600 { 00601 p(2) = v; 00602 break; 00603 } 00604 default: 00605 libmesh_error(); 00606 } 00607 } 00608 00609 // If we are in debug mode do a sanity check. Make sure 00610 // the point \p p on the reference element actually does 00611 // map to the point \p physical_point within a tolerance. 00612 #ifdef DEBUG 00613 /* 00614 const Point check = InfFE<Dim,T_radial,T_map>::map (inf_elem, p); 00615 const Point diff = physical_point - check; 00616 00617 if (diff.size() > tolerance) 00618 { 00619 libmesh_here(); 00620 libMesh::err << "WARNING: diff is " 00621 << diff.size() 00622 << std::endl; 00623 } 00624 */ 00625 #endif 00626 00627 00628 00629 00630 // Stop logging the map inversion. 00631 STOP_LOG("inverse_map()", "InfFE"); 00632 00633 return p; 00634 } 00635 00636 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00637 void InfFE<Dim,T_radial,T_map>::inverse_map (const Elem* elem, 00638 const std::vector<Point>& physical_points, 00639 std::vector<Point>& reference_points, 00640 const Real tolerance, 00641 const bool secure) 00642 { 00643 // The number of points to find the 00644 // inverse map of 00645 const std::size_t n_points = physical_points.size(); 00646 00647 // Resize the vector to hold the points 00648 // on the reference element 00649 reference_points.resize(n_points); 00650 00651 // Find the coordinates on the reference 00652 // element of each point in physical space 00653 for (unsigned int p=0; p<n_points; p++) 00654 reference_points[p] = 00655 InfFE<Dim,T_radial,T_map>::inverse_map (elem, physical_points[p], tolerance, secure, false); 00656 } 00657 00658 00659 00660 00661 //-------------------------------------------------------------- 00662 // Explicit instantiations using the macro from inf_fe_macro.h 00663 //INSTANTIATE_INF_FE(1,CARTESIAN); 00664 00665 //INSTANTIATE_INF_FE(2,CARTESIAN); 00666 00667 //INSTANTIATE_INF_FE(3,CARTESIAN); 00668 00669 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,Point,inverse_map(const Elem*,const Point&,const Real,const bool,const bool)); 00670 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,Point,inverse_map(const Elem*,const Point&,const Real,const bool,const bool)); 00671 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,Point,inverse_map(const Elem*,const Point&,const Real,const bool,const bool)); 00672 00673 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,inverse_map(const Elem*,const std::vector<Point>&,std::vector<Point>&,const Real, const bool)); 00674 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,inverse_map(const Elem*,const std::vector<Point>&,std::vector<Point>&,const Real, const bool)); 00675 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,inverse_map(const Elem*,const std::vector<Point>&,std::vector<Point>&,const Real, const bool)); 00676 00677 00678 } // namespace libMesh 00679 00680 00681 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00682
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: