$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 // C++ includes 00020 #include <algorithm> // for std::fill 00021 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00022 #include <cmath> // for sqrt 00023 00024 00025 // Local Includes 00026 #include "libmesh/libmesh_common.h" 00027 #include "libmesh/kelly_error_estimator.h" 00028 #include "libmesh/error_vector.h" 00029 #include "libmesh/fe_base.h" 00030 #include "libmesh/libmesh_logging.h" 00031 #include "libmesh/elem.h" 00032 #include "libmesh/system.h" 00033 00034 #include "libmesh/dense_vector.h" 00035 #include "libmesh/tensor_tools.h" 00036 00037 namespace libMesh 00038 { 00039 00040 00041 00042 void 00043 KellyErrorEstimator::initialize(const System& system, 00044 ErrorVector&, 00045 bool) 00046 { 00047 // Hang onto the system - we may need it for variable names later. 00048 my_system = &system; 00049 00050 // We'll need gradients and normal vectors for flux jump computation 00051 fe_fine->get_dphi(); 00052 fe_fine->get_normals(); 00053 fe_coarse->get_dphi(); 00054 } 00055 00056 00057 00058 void 00059 KellyErrorEstimator::internal_side_integration () 00060 { 00061 Real error = 1.e-30; 00062 unsigned int n_qp = fe_fine->n_quadrature_points(); 00063 unsigned int n_fine_dofs = Ufine.size(); 00064 unsigned int n_coarse_dofs = Ucoarse.size(); 00065 00066 std::vector<std::vector<RealGradient> > dphi_coarse = fe_coarse->get_dphi(); 00067 std::vector<std::vector<RealGradient> > dphi_fine = fe_fine->get_dphi(); 00068 std::vector<Point> face_normals = fe_fine->get_normals(); 00069 std::vector<Real> JxW_face = fe_fine->get_JxW(); 00070 00071 for (unsigned int qp=0; qp != n_qp; ++qp) 00072 { 00073 // Calculate solution gradients on fine and coarse elements 00074 // at this quadrature point 00075 Gradient grad_fine, grad_coarse; 00076 for (unsigned int i=0; i != n_coarse_dofs; ++i) 00077 grad_coarse.add_scaled (dphi_coarse[i][qp], Ucoarse(i)); 00078 00079 for (unsigned int i=0; i != n_fine_dofs; ++i) 00080 grad_fine.add_scaled (dphi_fine[i][qp], Ufine(i)); 00081 00082 // Find the jump in the normal derivative 00083 // at this quadrature point 00084 const Number jump = (grad_fine - grad_coarse)*face_normals[qp]; 00085 const Real jump2 = TensorTools::norm_sq(jump); 00086 00087 // Accumulate the jump integral 00088 error += JxW_face[qp] * jump2; 00089 } 00090 00091 // Add the h-weighted jump integral to each error term 00092 fine_error = 00093 error * fine_elem->hmax() * error_norm.weight(var); 00094 coarse_error = 00095 error * coarse_elem->hmax() * error_norm.weight(var); 00096 } 00097 00098 00099 bool 00100 KellyErrorEstimator::boundary_side_integration () 00101 { 00102 const std::string &var_name = my_system->variable_name(var); 00103 const unsigned int n_fine_dofs = Ufine.size(); 00104 00105 std::vector<std::vector<RealGradient> > dphi_fine = fe_fine->get_dphi(); 00106 std::vector<Point> face_normals = fe_fine->get_normals(); 00107 std::vector<Real> JxW_face = fe_fine->get_JxW(); 00108 std::vector<Point> qface_point = fe_fine->get_xyz(); 00109 00110 // The reinitialization also recomputes the locations of 00111 // the quadrature points on the side. By checking if the 00112 // first quadrature point on the side is on a flux boundary 00113 // for a particular variable, we will determine if the whole 00114 // element is on a flux boundary (assuming quadrature points 00115 // are strictly contained in the side). 00116 if (this->_bc_function(*my_system, qface_point[0], var_name).first) 00117 { 00118 const Real h = fine_elem->hmax(); 00119 00120 // The number of quadrature points 00121 const unsigned int n_qp = fe_fine->n_quadrature_points(); 00122 00123 // The error contribution from this face 00124 Real error = 1.e-30; 00125 00126 // loop over the integration points on the face. 00127 for (unsigned int qp=0; qp<n_qp; qp++) 00128 { 00129 // Value of the imposed flux BC at this quadrature point. 00130 const std::pair<bool,Real> flux_bc = 00131 this->_bc_function(*my_system, qface_point[qp], var_name); 00132 00133 // Be sure the BC function still thinks we're on the 00134 // flux boundary. 00135 libmesh_assert_equal_to (flux_bc.first, true); 00136 00137 // The solution gradient from each element 00138 Gradient grad_fine; 00139 00140 // Compute the solution gradient on element e 00141 for (unsigned int i=0; i != n_fine_dofs; i++) 00142 grad_fine.add_scaled (dphi_fine[i][qp], Ufine(i)); 00143 00144 // The difference between the desired BC and the approximate solution. 00145 const Number jump = flux_bc.second - grad_fine*face_normals[qp]; 00146 00147 // The flux jump squared. If using complex numbers, 00148 // TensorTools::norm_sq(z) returns |z|^2, where |z| is the modulus of z. 00149 const Real jump2 = TensorTools::norm_sq(jump); 00150 00151 // Integrate the error on the face. The error is 00152 // scaled by an additional power of h, where h is 00153 // the maximum side length for the element. This 00154 // arises in the definition of the indicator. 00155 error += JxW_face[qp]*jump2; 00156 00157 } // End quadrature point loop 00158 00159 fine_error = error*h*error_norm.weight(var); 00160 00161 return true; 00162 } // end if side on flux boundary 00163 return false; 00164 } 00165 00166 00167 00168 void 00169 KellyErrorEstimator::attach_flux_bc_function (std::pair<bool,Real> fptr(const System& system, 00170 const Point& p, 00171 const std::string& var_name)) 00172 { 00173 _bc_function = fptr; 00174 00175 // We may be turning boundary side integration on or off 00176 if (fptr) 00177 integrate_boundary_sides = true; 00178 else 00179 integrate_boundary_sides = false; 00180 } 00181 00182 } // namespace libMesh